API Reference

Classes

class AnnulusBounds : public Acts::DiscBounds

Class that implements a (potentially asymmetric) bounds with difference between surface bound center and surface coordinate center.

These bounds combine two different systems:

  • module system : radial bounds centred on the moduleOrigin

  • strip system : phi bounds centred on the stripOrigin

The measurement will be done in the strip system, with r/phi local coordinates.

Public Types

enum BoundValues

Values:

enumerator eMinR
enumerator eMaxR
enumerator eMinPhiRel
enumerator eMaxPhiRel
enumerator eAveragePhi
enumerator eOriginX
enumerator eOriginY
enumerator eSize

Public Functions

AnnulusBounds() = delete
AnnulusBounds(const AnnulusBounds &source) = default
AnnulusBounds(const std::array<double, eSize> &values) noexcept(false)

Constructor - from parameters array.

Parameters:

values – The parameter array

inline AnnulusBounds(double minR, double maxR, double minPhiRel, double maxPhiRel, const Vector2 &moduleOrigin = {0, 0}, double avgPhi = 0) noexcept(false)

Default constructor from parameters.

Note

For morigin you need to actually calculate the cartesian offset

Parameters:
  • minR – inner radius, in module system

  • maxR – outer radius, in module system

  • minPhiRel – right angular edge, in strip system, rel to avgPhi

  • maxPhiRel – left angular edge, in strip system, rel to avgPhi

  • moduleOrigin – The origin offset between the two systems.

  • avgPhi – (Optional) internal rotation of this bounds object’s local frame

inline virtual double binningValuePhi() const final

Return a reference radius for binning.

inline virtual double binningValueR() const final

Return a reference radius for binning.

std::vector<Vector2> corners() const

This method returns the four corners of the bounds in polar coordinates Starting from the upper right (max R, pos locX) and proceeding clock-wise i.e.

(max R; pos locX), (min R; pos locX), (min R; neg loc X), (max R: neg locX)

inline virtual bool coversFullAzimuth() const final

Returns true for full phi coverage.

inline double get(BoundValues bValue) const

Access to the bound values.

Parameters:

bValue – the class nested enum for the array access

virtual bool inside(const Vector2 &lposition, const BoundaryTolerance &boundaryTolerance) const final

Inside check for the bounds object driven by the boundary check directive Each Bounds has a method inside, which checks if a LocalPosition is inside the bounds Inside can be called without/with tolerances.

Parameters:
  • lposition – Local position (assumed to be in right surface frame)

  • boundaryTolerance – boundary check directive

Returns:

boolean indicator for the success of this operation

inline virtual bool insideRadialBounds(double R, double tolerance = 0.) const final

Checks if this is inside the radial coverage given the a tolerance.

Vector2 moduleOrigin() const

Returns moduleOrigin, but rotated out, so averagePhi is already considered.

The module origin needs to consider the rotation introduced by averagePhi

Returns:

The origin of the local frame

inline double phiMax() const

Returns the left angular edge of the module.

Returns:

The left side angle

inline double phiMin() const

Returns the right angular edge of the module.

Returns:

The right side angle

inline virtual double rMax() const final

This method returns outer radius.

inline virtual double rMin() const final

This method returns inner radius.

virtual std::ostream &toStream(std::ostream &sl) const final

Outstream operator.

Parameters:

sl – is the ostream to be dumped into

inline virtual SurfaceBounds::BoundsType type() const final

Return the bounds type - for persistency optimization.

Returns:

is a BoundsType enum

inline virtual std::vector<double> values() const final

Return the bound values as dynamically sized vector.

Returns:

this returns a copy of the internal values

virtual std::vector<Vector2> vertices(unsigned int quarterSegments = 2u) const override

This method returns the xy coordinates of the four corners of the bounds in module coordinates (in x/y), and if quarterSegments is bigger or equal to 0, the curved part of the segment is included and approximated by the corresponding number of segments.

Starting from the upper right (max R, pos locX) and proceeding clock-wise i.e. (max R; pos locX), (min R; pos locX), (min R; neg loc X), (max R: neg locX)

Parameters:

quarterSegments – the number of segments used to approximate a quarter of a circle

Returns:

vector for vertices in 2D

class AnyCharge

Charge and momentum interpretation for arbitrarily charged particles.

Only a charge magnitude identical to zero is interpreted as representing a neutral particle. This avoids ambiguities that might arise from using an approximate comparison with an arbitrary epsilon.

Public Functions

inline constexpr AnyCharge(float absQ) noexcept

Construct with the magnitude of the input charge.

inline constexpr AnyCharge(Neutral) noexcept
inline constexpr AnyCharge(SinglyCharged) noexcept
inline constexpr float absQ() const noexcept
inline constexpr float extractCharge(ActsScalar qOverP) const noexcept
inline constexpr ActsScalar extractMomentum(ActsScalar qOverP) const noexcept
inline constexpr ActsScalar qOverP(ActsScalar momentum, float signedQ) const noexcept

Friends

inline friend constexpr friend bool operator== (AnyCharge lhs, AnyCharge rhs) noexcept

Compare for equality.

template<int NComponents, int PolyDegree>
class AtlasBetheHeitlerApprox

This class approximates the Bethe-Heitler distribution as a gaussian mixture.

To enable an approximation for continuous input variables, the weights, means and variances are internally parametrized as a Nth order polynomial.

Public Types

using Data = std::array<PolyData, NComponents>

Public Functions

inline constexpr AtlasBetheHeitlerApprox(const Data &lowData, const Data &highData, bool lowTransform, bool highTransform, double lowLimit = 0.1, double highLimit = 0.2, bool clampToRange = false)

Construct the Bethe-Heitler approximation description with two parameterizations, one for lower ranges, one for higher ranges.

Is it assumed that the lower limit of the high-x/x0 data is equal to the upper limit of the low-x/x0 data.

Parameters:
  • lowData – data for the lower x/x0 range

  • highData – data for the higher x/x0 range

  • lowTransform – whether the low data need to be transformed

  • highTransform – whether the high data need to be transformed

  • lowLimit – the upper limit for the low data

  • highLimit – the upper limit for the high data

  • clampToRange – whether to clamp the input x/x0 to the allowed range

inline auto mixture(ActsScalar x) const

Generates the mixture from the polynomials and reweights them, so that the sum of all weights is 1.

Parameters:

x – pathlength in terms of the radiation length

inline constexpr auto numComponents() const

Returns the number of components the returned mixture will have.

inline constexpr bool validXOverX0(ActsScalar x) const

Checks if an input is valid for the parameterization.

Parameters:

x – pathlength in terms of the radiation length

Public Static Functions

static inline auto loadFromFiles(const std::string &low_parameters_path, const std::string &high_parameters_path, double lowLimit = 0.1, double highLimit = 0.2, bool clampToRange = false)

Loads a parameterization from a file according to the Atlas file description.

Parameters:
  • low_parameters_path – Path to the foo.par file that stores the parameterization for low x/x0

  • high_parameters_path – Path to the foo.par file that stores the parameterization for high x/x0

  • lowLimit – the upper limit for the low x/x0-data

  • highLimit – the upper limit for the high x/x0-data

  • clampToRange – forwarded to constructor

struct PolyData

Public Members

std::array<ActsScalar, PolyDegree + 1> meanCoeffs
std::array<ActsScalar, PolyDegree + 1> varCoeffs
std::array<ActsScalar, PolyDegree + 1> weightCoeffs
class AtlasStepper

the AtlasStepper implementation for the

Public Types

using BoundState = std::tuple<BoundTrackParameters, Jacobian, double>
using Covariance = BoundSquareMatrix
using CurvilinearState = std::tuple<CurvilinearTrackParameters, Jacobian, double>
using Jacobian = BoundMatrix

Public Functions

inline explicit AtlasStepper(const Config &config)
inline explicit AtlasStepper(std::shared_ptr<const MagneticFieldProvider> bField)
inline double absoluteMomentum(const State &state) const

Absolute momentum accessor.

Parameters:

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

inline 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.

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] Correction for non-linearity effect during transform from free to bound

Returns:

A bound state:

  • the parameters at the surface

  • the stepwise jacobian towards it

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

Create and return a curvilinear state at the current position.

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

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

inline Vector3 direction(const State &state) const
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, then it takes the field from the cell.

Parameters:
  • state[inout] is the stepper 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.

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

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

inline State makeState(std::reference_wrapper<const GeometryContext> gctx, std::reference_wrapper<const MagneticFieldContext> mctx, const BoundTrackParameters &par, double ssize = std::numeric_limits<double>::max(), double stolerance = s_onSurfaceTolerance) const
inline Vector3 momentum(const State &state) const
inline std::string outputStepSize(const State &state) const

Output the Step Size - single component.

Parameters:

state[inout] The stepping state (thread-local cache)

inline double overstepLimit(const State&) const

Overstep limit.

inline const ParticleHypothesis &particleHypothesis(const State &state) const

Particle hypothesis.

Parameters:

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

inline Vector3 position(const State &state) const
template<typename propagator_state_t, typename navigator_t>
inline bool prepareCurvilinearState([[maybe_unused]] propagator_state_t &prop_state, [[maybe_unused]] const navigator_t &navigator) const

If necessary fill additional members needed for curvilinearState.

Compute path length derivatives in case they have not been computed yet, which is the case if no step has been executed yet.

Parameters:
  • prop_state[inout] State that will be presented as BoundState

  • navigator[in] the navigator of the propagation

Returns:

true if nothing is missing after this call, false otherwise.

inline double qOverP(const State &state) const
inline void releaseStepSize(State &state, ConstrainedStep::Type stype) const

Release the Step size.

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

  • stype[in] The step size type to be released

inline void resetState(State &state, const BoundVector &boundParams, const BoundSquareMatrix &cov, const Surface &surface, 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] Reset state will be on this surface

  • stepSize[in] Step size

inline void setIdentityJacobian(State &state) const

Method that reset the Jacobian to the Identity for when no bound state are available.

Parameters:

state[inout] State of the stepper

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

Perform the actual step on the state.

Parameters:

state – is the provided stepper state (caller keeps thread locality)

inline double time(const State &state) const

Time access.

inline void transportCovarianceToBound(State &state, const Surface &surface, const 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.

Parameters:
  • state[inout] State of the stepper

  • surface[in] is the surface to which the covariance is forwarded 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

inline void update(State &state, const FreeVector &parameters, const BoundVector &boundParams, const Covariance &covariance, const Surface &surface) const

The state update method.

Parameters:
  • state[inout] The stepper state for

  • parameters[in] The new free track parameters at start

  • boundParams[in] Corresponding bound parameters

  • covariance[in] The updated covariance matrix

  • surface[in] The surface used to update the pVector

inline void update(State &state, const Vector3 &uposition, const Vector3 &udirection, double qop, double time) const

Method to update momentum, direction and p.

Parameters:
  • state – The state object

  • uposition – the updated position

  • udirection – the updated direction

  • qop – the updated momentum value

  • time – the update time

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

Update step size.

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

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

Update step size - explicitly with a double.

Parameters:
  • state[inout] 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?

inline Intersection3D::Status updateSurfaceStatus(State &state, const Surface &surface, std::uint8_t index, Direction navDir, const BoundaryTolerance &boundaryTolerance, ActsScalar surfaceTolerance = s_onSurfaceTolerance, const Logger &logger = getDummyLogger()) const

Update surface status.

This method intersect 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[inout] The stepping state (thread-local cache)

  • surface[in] The surface provided

  • index[in] The surface intersection index

  • navDir[in] The navigation direction

  • boundaryTolerance[in] The boundary check for this status update

  • surfaceTolerance[in] Surface tolerance used for intersection

  • logger[in] Logger instance to use

struct Config

Public Members

std::shared_ptr<const MagneticFieldProvider> bField
struct Options : public Acts::StepperPlainOptions

Public Functions

inline void setPlainOptions(const StepperPlainOptions &options)
struct State

Nested State struct for the local caching.

Public Functions

State() = delete

Default constructor - deleted.

template<typename Parameters>
inline State(const GeometryContext &gctx, MagneticFieldProvider::Cache fieldCacheIn, const Parameters &pars, double ssize = std::numeric_limits<double>::max(), double stolerance = s_onSurfaceTolerance)

Constructor.

Template Parameters:

Type – of TrackParameters

Parameters:
  • gctx[in] The geometry context tof this call

  • fieldCacheIn[in] The magnetic field cache for this call

  • pars[in] Input parameters

  • ssize[in] the steps size limitation

  • stolerance[in] is the stepping tolerance

Public Members

Covariance cov = Covariance::Zero()
const Covariance *covariance = nullptr
bool covTransport = false
bool debug = false

Debug output the string where debug messages are stored (optionally)

std::size_t debugMsgWidth = 50
std::size_t debugPfxWidth = 30

buffer & formatting for consistent output

std::string debugString = ""
Vector3 field
MagneticFieldProvider::Cache fieldCache

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

See step() code for details.

std::reference_wrapper<const GeometryContext> geoContext

Cache the geometry context.

double jacobian[eBoundSize * eBoundSize] = {}
double maxPathLength = 0
bool mcondition = false
bool needgradient = false
bool newfield = true
std::size_t nSteps = 0

Total number of performed steps.

std::size_t nStepTrials = 0

Totoal number of attempted steps.

double parameters[eBoundSize] = {0., 0., 0., 0., 0., 0.}

Storage pattern of pVector /dL0 /dL1 /dPhi /the /dCM /dT X ->P[0] dX / P[ 8] P[16] P[24] P[32] P[40] P[48] Y ->P[1] dY / P[ 9] P[17] P[25] P[33] P[41] P[49] Z ->P[2] dZ / P[10] P[18] P[26] P[34] P[42] P[50] T ->P[3] dT/ P[11] P[19] P[27] P[35] P[43] P[51] Ax ->P[4] dAx/ P[12] P[20] P[28] P[36] P[44] P[52] Ay ->P[5] dAy/ P[13] P[21] P[29] P[37] P[45] P[53] Az ->P[6] dAz/ P[14] P[22] P[30] P[38] P[46] P[54] CM ->P[7] dCM/ P[15] P[23] P[31] P[39] P[47] P[55] Cache: P[56] - P[59].

ParticleHypothesis particleHypothesis
double pathAccumulated = 0.
double previousStepSize = 0.
std::array<double, 60> pVector = {}
bool state_ready = false
double step = 0
ConstrainedStep stepSize
double tolerance = s_onSurfaceTolerance

The tolerance for the stepping.

bool useJacobian = false
template<AxisType type, AxisBoundaryType bdt = AxisBoundaryType::Open>
class Axis

calculate bin indices from a given binning structure

This class provides some basic functionality for calculating bin indices for a given binning configuration. Both equidistant as well as variable binning structures are supported.

Bin intervals are defined such that the lower bound is closed and the upper bound is open.

Template Parameters:

equidistant – flag whether binning is equidistant (true) or not (false)

class BinUtility

The BinUtility class that translated global and local position into a bins of a BinnedArray, most performant is equidistant binning without a transform, however, optionally a transform can be provided, e.g.

for binning on shifted object, the transform is usually shared with the geometric object the Array is defined on, for performance reasons, also the inverse transform is stored.

Public Functions

inline BinUtility()

Constructor for equidistant.

BinUtility(BinUtility &&sbu) = default
inline BinUtility(const BinningData &bData, const Transform3 &tForm = Transform3::Identity())

Constructor from BinningData directly.

Parameters:
  • bData – is the provided binning data

  • tForm – is the (optional) transform

BinUtility(const BinUtility &sbu) = default

Copy constructor.

Parameters:

sbu – is the source bin utility

inline explicit BinUtility(const Transform3 &tForm)

Constructor with only a Transform3.

Parameters:

tForm – is the local to global transform

inline BinUtility(std::size_t bins, float min, float max, BinningOption opt = open, BinningValue value = BinningValue::binX, const Transform3 &tForm = Transform3::Identity())

Constructor for equidistant.

Parameters:
  • bins – is the number of bins

  • min – in the minimal value

  • max – is the maximal value

  • opt – is the binning option : open, closed

  • value – is the binninb value : binX, binY, binZ, etc.

  • tForm – is the (optional) transform

inline BinUtility(std::vector<float> &bValues, BinningOption opt = open, BinningValue value = BinningValue::binPhi, const Transform3 &tForm = Transform3::Identity())

Constructor for arbitrary.

Parameters:
  • bValues – is the boundary values of the binning

  • opt – is the binning option : open, closed

  • value – is the binninb value : binX, binY, binZ, etc.

  • tForm – is the (optional) transform

~BinUtility() = default

Virtual Destructor.

inline std::size_t bin(const Vector2 &lposition, std::size_t ba = 0) const

Bin from a 2D vector (following local parameters defintitions)

  • no optional transform applied

  • USE WITH CARE !!

You need to make sure that the local position is actually in the binning frame of the BinUtility

Parameters:
  • lposition – is the local position to be set

  • ba – is the bin dimension

Returns:

bin calculated from local

inline std::size_t bin(const Vector3 &position, std::size_t ba = 0) const

Bin from a 3D vector (already in binning frame)

Parameters:
  • position – is the 3D position to be evaluated

  • ba – is the bin dimension

Returns:

is the bin value

inline const std::vector<BinningData> &binningData() const

Return the binning data vector.

inline BinningValue binningValue(std::size_t ba = 0) const

The type/value of the binning.

Parameters:

ba – is the binaccessor

Returns:

the binning value of the accessor entry

inline std::size_t bins() const

Return the total number of bins.

inline std::size_t bins(std::size_t ba) const

Number of bins.

Parameters:

ba – is the binaccessor

Returns:

std::size_t is the bins of the accessor entry

inline std::array<std::size_t, 3> binTriple(const Vector3 &position) const

Bin-triple fast access.

  • calculate the bin triple with one transform

Parameters:

position – is the 3D position to be evaluated

Returns:

is the bin value in 3D

inline std::size_t dimensions() const

First bin maximal value.

Returns:

the dimension of the binning data

inline bool inside(const Vector3 &position) const

Check if bin is inside from Vector2 - optional transform applied.

Parameters:

position – is the global position to be evaluated

Returns:

is a boolean check

inline std::size_t max(std::size_t ba = 0) const

First bin maximal value.

Parameters:

ba – is the binaccessor

Returns:

std::size_t is the maximal bin of the accessor entry

inline int nextDirection(const Vector3 &position, const Vector3 &direction, std::size_t ba = 0) const

Return the other direction for fast interlinking.

Parameters:
  • position – is the global position for the next search

  • direction – is the global position for the next search

  • ba – is the bin accessor

Returns:

the next bin

inline BinUtility &operator+=(const BinUtility &gbu)

Operator+= to make multidimensional BinUtility.

Parameters:

gbu – is the additional BinUtility to be chosen

BinUtility &operator=(BinUtility&&) = default
inline BinUtility &operator=(const BinUtility &sbu)

Assignment operator.

Parameters:

sbu – is the source bin utility

inline bool operator==(const BinUtility &other) const

Equality operator.

inline std::size_t serialize(const std::array<std::size_t, 3> &bin) const

Serialize the bin triple.

  • this creates a simple std::size_t from a triple object

Parameters:

bin – is the bin to be serialized

inline std::ostream &toStream(std::ostream &sl, const std::string &indent = "") const

Output Method for std::ostream, to be overloaded by child classes.

Parameters:
  • sl – is the ostream to be dumped into

  • indent – the current indentation

Returns:

the input stream

inline std::string toString(const std::string &indent = "") const

Output into a string.

Parameters:

indent – the current indentation

Returns:

a string with the stream information

inline const Transform3 &transform() const

Transform applied to global positions before lookup.

Returns:

Shared pointer to transform

Friends

inline friend std::ostream &operator<<(std::ostream &sl, const BinUtility &bgen)

Overload of << operator for std::ostream for debug output.

class BinnedSurfaceMaterial : public Acts::ISurfaceMaterial

It extends the ISurfaceMaterial base class and is an array pf MaterialSlab.

This is not memory optimised as every bin holds one material property object.

Public Functions

BinnedSurfaceMaterial() = delete

Default Constructor - deleted.

BinnedSurfaceMaterial(BinnedSurfaceMaterial &&bsm) = default

Copy Move Constructor.

Parameters:

bsm – is the source object to be copied

BinnedSurfaceMaterial(const BinnedSurfaceMaterial &bsm) = default

Copy Constructor.

Parameters:

bsm – is the source object to be copied

BinnedSurfaceMaterial(const BinUtility &binUtility, MaterialSlabMatrix fullProperties, double splitFactor = 0., MappingType mappingType = MappingType::Default)

Explicit constructor with only full MaterialSlab, for two-dimensional binning.

The split factors:

  • 1. : oppositePre

  • 0. : alongPre ===> 1 Dimensional array

Parameters:
  • binUtility – defines the binning structure on the surface (copied)

  • fullProperties – is the vector of properties as recorded (moved)

  • splitFactor – is the pre/post splitting directive

  • mappingType – is the type of surface mapping associated to the surface

BinnedSurfaceMaterial(const BinUtility &binUtility, MaterialSlabVector fullProperties, double splitFactor = 0., MappingType mappingType = MappingType::Default)

Explicit constructor with only full MaterialSlab, for one-dimensional binning.

The split factors:

  • 1. : oppositePre

  • 0. : alongPre ===> 1 Dimensional array

Parameters:
  • binUtility – defines the binning structure on the surface (copied)

  • fullProperties – is the vector of properties as recorded (moved)

  • splitFactor – is the pre/post splitting directive

  • mappingType – is the type of surface mapping associated to the surface

~BinnedSurfaceMaterial() override = default

Destructor.

inline const BinUtility &binUtility() const

Return the BinUtility.

inline const MaterialSlabMatrix &fullMaterial() const

Retrieve the entire material slab matrix.

virtual const MaterialSlab &materialSlab(const Vector2 &lp) const final

Return method for full material description of the Surface.

  • from local coordinate on the surface

Parameters:

lp – is the local position used for the (eventual) lookup

Returns:

const MaterialSlab

virtual const MaterialSlab &materialSlab(const Vector3 &gp) const final

Return method for full material description of the Surface.

  • from the global coordinates

Parameters:

gp – is the global position used for the (eventual) lookup

Returns:

const MaterialSlab

virtual BinnedSurfaceMaterial &operator*=(double scale) final

Scale operator.

Parameters:

scale – is the scale factor for the full material

BinnedSurfaceMaterial &operator=(BinnedSurfaceMaterial &&bsm) = default

Assignment Move operator.

BinnedSurfaceMaterial &operator=(const BinnedSurfaceMaterial &bsm) = default

Assignment operator.

virtual std::ostream &toStream(std::ostream &sl) const final

Output Method for std::ostream, to be overloaded by child classes.

class BinningData

This class holds all the data necessary for the bin calculation.

phi has a very particular behaviour:

  • there’s the change around +/- PI

  • it can be multiplicative or additive multiplicative : each major bin has the same sub structure i.e. first binnning

structure is equidistant additive : sub structure replaces one bin (and one bin only)

Public Functions

BinningData() = default
inline BinningData(BinningOption bOption, BinningValue bValue, const std::vector<float> &bBoundaries, std::unique_ptr<const BinningData> sBinData = nullptr)

Constructor for non-equidistant binning.

Parameters:
  • bOption – is the binning option : open / closed

  • bValue – is the binning value : binX, binY, etc.

  • bBoundaries – are the bin boundaries

  • sBinData – is (optional) sub structure

inline BinningData(BinningOption bOption, BinningValue bValue, std::size_t bBins, float bMin, float bMax, std::unique_ptr<const BinningData> sBinData = nullptr, bool sBinAdditive = false)

Constructor for equidistant binning and optional sub structure can be multiplicative or additive.

Parameters:
  • bOption – is the binning option : open, closed

  • bValue – is the binning value: binX, binY, etc.

  • bBins – is number of equidistant bins

  • bMin – is the minimum value

  • bMax – is the maximum value

  • sBinData – is (optional) sub structure

  • sBinAdditive – is the prescription for the sub structure

inline BinningData(BinningValue bValue, float bMin, float bMax)

Constructor for 0D binning.

Parameters:
  • bValue – is the binning value: binX, binY, etc.

  • bMin – is the minimum value

  • bMax – is the maximum value

inline BinningData(const BinningData &bdata)

Copy constructor.

Parameters:

bdata – is the source object

~BinningData() = default
inline std::size_t bins() const

Return the number of bins - including sub bins.

inline const std::vector<float> &boundaries() const

Return the boundaries - including sub boundaries.

Returns:

vector of floats indicating the boundary values

inline float center(std::size_t bin) const

Get the center value of a bin.

Parameters:

bin – is the bin for which the center value is requested

Returns:

float value according to the bin center

inline float centerValue(std::size_t bin) const

access to the center value this uses the bin boundary vector, it also works with sub structure

Parameters:

bin – is the bin for which the value is requested, if bin > nbins it is set to max

Returns:

the center value of the bin is given

inline bool inside(const Vector2 &lposition) const

Check if bin is inside from Vector2.

Parameters:

lposition – is the search position in global coordinated

Returns:

boolean if this is inside() method is true

inline bool inside(const Vector3 &position) const

Check if bin is inside from Vector3.

Parameters:

position – is the search position in global coordinated

Returns:

boolean if this is inside() method is true

inline int nextDirection(const Vector3 &position, const Vector3 &dir) const

Layer next direction is needed.

Parameters:
  • position – is the start search position

  • dir – is the direction

Returns:

integer that indicates which direction to move

inline BinningData &operator=(const BinningData &bdata)

Assignment operator.

Parameters:

bdata – is the source object

inline bool operator==(const BinningData &bData) const

Equality operator.

Parameters:

bData – is the binning data to be checked against

Returns:

a boolean indicating if they are the same

inline std::size_t search(float value) const

Generic search - forwards to correct function pointer.

Parameters:

value – is the searchvalue as float

Returns:

bin according tot this

inline std::size_t searchGlobal(const Vector3 &position) const

Generic search from a 3D position &#8212; corresponds to global coordinate schema.

Parameters:

position – is the search position in global coordinated

Returns:

bin according tot this

inline std::size_t searchLocal(const Vector2 &lposition) const

Generic search from a 2D position &#8212; corresponds to local coordinate schema.

Parameters:

lposition – is the search position in local coordinated

Returns:

bin according tot this

inline std::size_t searchWithSubStructure(float value) const

Generic search with sub structure.

  • forwards to correct function pointer

Parameters:

value – is the searchvalue as float

Returns:

bin according tot this

inline std::string toString(const std::string &indent = "") const

String screen output method.

Parameters:

indent – the current indentation

Returns:

a string containing the screen information

inline float value(const Vector2 &lposition) const

Take the right float value.

Parameters:

lposition – assumes the correct local position expression

Returns:

float value according to the binning setup

inline float value(const Vector3 &position) const

Take the right float value.

Parameters:

position – is the global position

Returns:

float value according to the binning setup

inline float width(std::size_t bin) const

Get the width of a bin.

Parameters:

bin – is the bin for which the width is requested

Returns:

float value of width

Public Members

BinningValue binvalue = {}

binning value: binX, binY, binZ, binR …

float max = {}

maximum value

float min = {}

minimum value

BinningOption option = {}

binning option: open, closed

float step = {}

binning step

bool subBinningAdditive = {}

sub structure: additive or multipicative

std::unique_ptr<const BinningData> subBinningData

sub structure: describe some sub binning

BinningType type = {}

binning type: equidistant, arbitrary

bool zdim = {}

zero dimensional binning : direct access

class BoundaryTolerance

Variant-like type to capture different types of boundary tolerances.

Since our track hypothesis comes with uncertainties, we sometimes need to check if the track is not just within the boundary of the surface but also within a certain tolerance. This class captures different parameterizations of such tolerances. The surface class will then use these tolerances to check if a ray is within the boundary+tolerance of the surface.

Different types of boundary tolerances implemented:

  • Infinite: Infinite tolerance i.e. no boundary check will be performed.

  • None: No tolerance i.e. exact boundary check will be performed.

  • AbsoluteBound: Absolute tolerance in bound coordinates. The tolerance is defined as a pair of absolute values for the bound coordinates. Only if both coordinates are within the tolerance, the boundary check is considered as passed.

  • AbsoluteCartesian: Absolute tolerance in Cartesian coordinates. The tolerance is defined as a pair of absolute values for the Cartesian coordinates. The transformation to Cartesian coordinates can be done via the Jacobian for small distances. Only if both coordinates are within the tolerance, the boundary check is considered as passed.

  • AbsoluteEuclidean: Absolute tolerance in Euclidean distance. The tolerance is defined as a single absolute value for the Euclidean distance. The Euclidean distance can be calculated via the local bound Jacobian and the bound coordinate residual. If the distance is within the tolerance, the boundary check is considered as passed.

  • Chi2Bound: Chi2 tolerance in bound coordinates. The tolerance is defined as a maximum chi2 value and a weight matrix, which is the inverse of the bound covariance matrix. The chi2 value is calculated from the bound coordinates residual and the weight matrix. If the chi2 value is below the maximum chi2 value, the boundary check is considered as passed.

The bound coordinates residual is defined as the difference between the point checked and the closest point on the boundary. The Jacobian is the derivative of the bound coordinates with respect to the Cartesian coordinates.

Public Types

using Variant = std::variant<Infinite, None, AbsoluteBound, AbsoluteCartesian, AbsoluteEuclidean, Chi2Bound>

Underlying variant type.

Public Functions

BoundaryTolerance(const AbsoluteBound &AbsoluteBound)

Construct with absolute tolerance in bound coordinates.

BoundaryTolerance(const AbsoluteCartesian &absoluteCartesian)

Construct with absolute tolerance in Cartesian coordinates.

BoundaryTolerance(const AbsoluteEuclidean &absoluteEuclidean)

Construct with absolute tolerance in Euclidean distance.

BoundaryTolerance(const Chi2Bound &Chi2Bound)

Construct with chi2 tolerance in bound coordinates.

BoundaryTolerance(const Infinite &infinite)

Construct with infinite tolerance.

BoundaryTolerance(const None &none)

Construct with no tolerance.

BoundaryTolerance(Variant variant)

Construct from variant.

AbsoluteBound asAbsoluteBound(bool isCartesian = false) const

Get the tolerance as absolute bound.

std::optional<AbsoluteBound> asAbsoluteBoundOpt(bool isCartesian = false) const

Get the tolerance as absolute bound if possible.

const AbsoluteCartesian &asAbsoluteCartesian() const

Get the tolerance as absolute Cartesian.

const AbsoluteEuclidean &asAbsoluteEuclidean() const

Get the tolerance as absolute Euclidean.

const Chi2Bound &asChi2Bound() const

Get the tolerance as chi2 bound.

SquareMatrix2 getMetric(const std::optional<SquareMatrix2> &jacobian) const

Get the metric for the tolerance.

bool hasAbsoluteBound(bool isCartesian = false) const

Check if the tolerance is absolute with bound coordinates.

bool hasAbsoluteCartesian() const

Check if the tolerance is absolute with Cartesian coordinates.

bool hasAbsoluteEuclidean() const

Check if the tolerance is absolute with Euclidean distance.

bool hasChi2Bound() const

Check if the tolerance is chi2 with bound coordinates.

bool hasMetric(bool hasJacobian) const

Check if there is a metric assigned with this tolerance.

bool hasTolerance() const

Check if any tolerance is set.

bool isInfinite() const

Check if the tolerance is infinite.

bool isNone() const

Check if the is no tolerance.

bool isTolerated(const Vector2 &distance, const std::optional<SquareMatrix2> &jacobianOpt) const

Check if the distance is tolerated.

struct AbsoluteBound

Absolute tolerance in bound coordinates.

Public Functions

AbsoluteBound() = default
inline AbsoluteBound(double tolerance0_, double tolerance1_)

Public Members

double tolerance0 = {}
double tolerance1 = {}
struct AbsoluteCartesian

Absolute tolerance in Cartesian coordinates.

Public Functions

AbsoluteCartesian() = default
inline AbsoluteCartesian(double tolerance0_, double tolerance1_)

Public Members

double tolerance0 = {}
double tolerance1 = {}
struct AbsoluteEuclidean

Absolute tolerance in Euclidean distance.

Public Functions

AbsoluteEuclidean() = default
inline explicit AbsoluteEuclidean(double tolerance_)

Public Members

double tolerance = {}
struct Chi2Bound

Chi2 tolerance in bound coordinates.

Public Functions

Chi2Bound() = default
inline Chi2Bound(const SquareMatrix2 &weight_, double maxChi2_)

Public Members

double maxChi2 = {}
SquareMatrix2 weight = SquareMatrix2::Identity()
struct Infinite

Infinite tolerance i.e. no boundary check.

struct None

No tolerance i.e. exact boundary check.

class CalibrationContext : public ContextType

This is the central definition of the Acts payload object regarding detector calibration.

It is propagated through the code to allow for event/thread dependent calibration

template<typename propagator_t, typename track_container_t>
class CombinatorialKalmanFilter

Combinatorial Kalman filter to find tracks.

The CombinatorialKalmanFilter contains an Actor and a Sequencer sub-class. The Sequencer has to be part of the Navigator of the Propagator in order to initialize and provide the measurement surfaces.

The Actor is part of the Propagation call and does the Kalman update. Updater and Calibrator are given to the Actor for further use:

  • The Updater is the implemented kalman updater formalism, it runs via a visitor pattern through the measurements.

Measurements are not required to be ordered for the CombinatorialKalmanFilter, measurement ordering needs to be figured out by the navigation of the propagator.

The void components are provided mainly for unit testing.

Template Parameters:

propagator_t – Type of the propagator

Public Functions

CombinatorialKalmanFilter() = delete

Default constructor is deleted.

inline CombinatorialKalmanFilter(propagator_t pPropagator, std::unique_ptr<const Logger> _logger = getDefaultLogger("CKF", Logging::INFO))

Constructor from arguments.

Combinatorial Kalman Filter implementation, calls the Kalman filter.

Note

The input measurements are given in the form of SourceLinks. It’s calibrator_t's job to turn them into calibrated measurements used in the track finding.

Template Parameters:
  • source_link_iterator_t – Type of the source link iterator

  • start_parameters_t – Type of the initial parameters

  • parameters_t – Type of parameters used for local parameters

Parameters:
  • initialParameters – The initial track parameters

  • tfOptions – CombinatorialKalmanFilterOptions steering the track finding

  • trackContainer – Track container in which to store the results

Returns:

a container of track finding result for all the initial track parameters

Combinatorial Kalman Filter implementation, calls the Kalman filter.

Note

The input measurements are given in the form of SourceLinks. It’s calibrator_t's job to turn them into calibrated measurements used in the track finding.

Template Parameters:
  • source_link_iterator_t – Type of the source link iterator

  • start_parameters_t – Type of the initial parameters

  • parameters_t – Type of parameters used for local parameters

Parameters:
  • initialParameters – The initial track parameters

  • tfOptions – CombinatorialKalmanFilterOptions steering the track finding

  • trackContainer – Track container in which to store the results

  • rootBranch – The track to be used as the root branch

Returns:

a container of track finding result for all the initial track parameters

class ConeBounds : public Acts::SurfaceBounds

Bounds for a conical surface, the opening angle is stored in \( \tan(\alpha) \) and always positively defined.

The cone can open to both sides, steered by \( z_min \) and \( z_max \).

../_images/ConeBounds.gif

Public Types

enum BoundValues

Values:

enumerator eAlpha
enumerator eMinZ
enumerator eMaxZ
enumerator eHalfPhiSector
enumerator eAveragePhi
enumerator eSize

Public Functions

ConeBounds() = delete
ConeBounds(const std::array<double, eSize> &values) noexcept(false)

Constructor - from parameters array.

Parameters:

values – The parameter array

ConeBounds(double alpha, bool symm, double halfphi = M_PI, double avphi = 0.) noexcept(false)

Constructor - open cone with alpha, by default a full cone but optionally can make a conical section.

Parameters:
  • alpha – is the opening angle of the cone

  • symm – is the boolean indicating if the cone is symmetric in +/- z

  • halfphi – is the half opening angle (default is pi)

  • avphi – is the phi value around which the bounds are opened (default=0)

ConeBounds(double alpha, double minz, double maxz, double halfphi = M_PI, double avphi = 0.) noexcept(false)

Constructor - open cone with alpha, minz and maxz, by default a full cone but can optionally make it a conical section.

Parameters:
  • alpha – is the opening angle of the cone

  • minz – cone expanding from minimal z

  • maxz – cone expanding to maximal z

  • halfphi – is the half opening angle (default is pi)

  • avphi – is the phi value around which the bounds are opened (default=0)

~ConeBounds() override = default
inline double get(BoundValues bValue) const

Access to the bound values.

Parameters:

bValue – the class nested enum for the array access

virtual bool inside(const Vector2 &lposition, const BoundaryTolerance &boundaryTolerance = BoundaryTolerance::None()) const final

inside method for local position

Parameters:
  • lposition – is the local position to be checked

  • boundaryTolerance – is the boundary check directive

Returns:

is a boolean indicating if the position is inside

inline double r(double z) const

Return the radius at a specific z values.

Parameters:

z – is the z value for which r is requested

Returns:

is the r value associated with z

inline double tanAlpha() const

Return tangent of alpha (pre-computed)

virtual std::ostream &toStream(std::ostream &sl) const final

Output Method for std::ostream.

Parameters:

sl – is the ostrea into which the dump is done

Returns:

is the input object

virtual BoundsType type() const final

Return the bounds type - for persistency optimization.

Returns:

is a BoundsType enum

inline virtual std::vector<double> values() const final

Return the bound values as dynamically sized vector.

Returns:

this returns a copy of the internal values

class ConeLayer : public virtual Acts::ConeSurface, public Acts::Layer

Class to describe a conical detector layer for tracking, it inhertis from both, Layer base class and ConeSurface class.

Public Functions

ConeLayer() = delete
ConeLayer(const ConeLayer &cla) = delete
~ConeLayer() override = default
ConeLayer &operator=(const ConeLayer&) = delete
virtual const ConeSurface &surfaceRepresentation() const override

Transforms the layer into a Surface representation for extrapolation.

virtual ConeSurface &surfaceRepresentation() override

Public Static Functions

static inline MutableLayerPtr create(const Transform3 &transform, std::shared_ptr<const ConeBounds> cbounds, std::unique_ptr<SurfaceArray> surfaceArray, double thickness = 0., std::unique_ptr<ApproachDescriptor> ad = nullptr, LayerType laytyp = Acts::active)

Factory for shared layer.

Parameters:
  • transform – is the 3D transform that positions the layer in 3D frame

  • cbounds – is the conical bound description

  • surfaceArray – is the array of sensitive surfaces

  • thickness – is the layer thickness along the normal axis

  • ad – is the approach descriptor for navigation towards the layer

  • laytyp – is the layer type

Returns:

is a shared pointer to a layer

class ConeSurface : public Acts::RegularSurface

Class for a conical surface in the Tracking geometry.

It inherits from Surface.

The ConeSurface is special since no corresponding Track parameters exist since they’re numerical instable at the tip of the cone. Propagations to a cone surface will be returned in curvilinear coordinates.

Subclassed by Acts::ConeLayer

Public Functions

ConeSurface() = delete
~ConeSurface() override = default
virtual AlignmentToPathMatrix alignmentToPathDerivative(const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction) const final

Calculate the derivative of path length at the geometry constraint or point-of-closest-approach w.r.t.

alignment parameters of the surface (i.e. local frame origin in global 3D Cartesian coordinates and its rotation represented with extrinsic Euler angles)

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • position – global 3D position

  • direction – global 3D momentum direction

Returns:

Derivative of path length w.r.t. the alignment parameters

virtual Vector3 binningPosition(const GeometryContext &gctx, BinningValue bValue) const final

The binning position method - is overloaded for r-type binning.

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • bValue – defines the type of binning applied in the global frame

Returns:

The return type is a vector for positioning in the global frame

virtual const ConeBounds &bounds() const final

This method returns the ConeBounds by reference.

Result<Vector2> globalToLocal(const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, double tolerance = s_onSurfaceTolerance) const final

Convert a global position to a local one this is the most generic interface, which is implemented by all surfaces.

Note

The position is required to be on-surface, which is indicated by the Result return value.

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • position – is the global position to be converted

  • direction – is the direction of the local position (ignored for RegularSurface)

  • tolerance – is the tolerance for the on-surface check

Returns:

Result type containing local position by value

virtual Result<Vector2> globalToLocal(const GeometryContext &gctx, const Vector3 &position, double tolerance = s_onSurfaceTolerance) const final

Global to local transformation.

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • position – is the global position to be transformed

  • tolerance – optional tolerance within which a point is considered valid on surface

Returns:

a Result<Vector2> which can be !ok() if the operation fails

Result<Vector2> globalToLocal(const GeometryContext &gctx, const Vector3 &position, double tolerance = s_onSurfaceTolerance) const = 0

Convert a global position to a local one.

Note

The position is required to be on-surface, which is indicated by the Result return value.

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • position – is the global position to be converted

  • tolerance – is the tolerance for the on-surface check

Returns:

Result type containing local position by value

virtual SurfaceMultiIntersection intersect(const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, const BoundaryTolerance &boundaryTolerance = BoundaryTolerance::Infinite(), double tolerance = s_onSurfaceTolerance) const final

Straight line intersection schema from position/direction.

If possible returns both solutions for the cylinder

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • position – The position to start from

  • direction – The direction at start

  • boundaryTolerance – the Boundary Check

  • tolerance – the tolerance used for the intersection

Returns:

SurfaceMultiIntersection object (contains intersection & surface)

virtual ActsMatrix<2, 3> localCartesianToBoundLocalDerivative(const GeometryContext &gctx, const Vector3 &position) const final

Calculate the derivative of bound track parameters local position w.r.t.

position in local 3D Cartesian coordinates

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • position – The position of the parameters in global

Returns:

Derivative of bound local position w.r.t. position in local 3D cartesian coordinates

virtual Vector3 localToGlobal(const GeometryContext &gctx, const Vector2 &lposition) const final

Local to global transformation.

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • lposition – is the local position to be transformed

Returns:

The global position by value

Vector3 localToGlobal(const GeometryContext &gctx, const Vector2 &lposition) const = 0

Local to global transformation.

Parameters:
  • gctx – The current geometry context object, e.g. alignment

  • lposition – local 2D position in specialized surface frame

Returns:

The global position by value

Vector3 localToGlobal(const GeometryContext &gctx, const Vector2 &lposition, const Vector3 &