Namespace Acts
-
namespace Acts
Note
This file is foreseen for the
Geometry
module to replaceExtent
Typedefs
-
using ActsDynamicMatrix = Eigen::Matrix<ActsScalar, Eigen::Dynamic, Eigen::Dynamic>
-
using ActsDynamicVector = Eigen::Matrix<ActsScalar, Eigen::Dynamic, 1>
-
template<unsigned int kRows, unsigned int kCols>
using ActsMatrix = Eigen::Matrix<ActsScalar, kRows, kCols>
-
using ActsScalar = double
Common scalar (floating point type used for the default algebra types.
Defaults to
double
but can be customized by the user.
-
template<unsigned int kSize>
using ActsSquareMatrix = Eigen::Matrix<ActsScalar, kSize, kSize>
-
template<unsigned int kSize>
using ActsVector = Eigen::Matrix<ActsScalar, kSize, 1>
-
using AlignmentMatrix = ActsMatrix<eAlignmentSize, eAlignmentSize>
-
using AlignmentRowVector = ActsMatrix<1, eAlignmentSize>
-
using AlignmentToBoundMatrix = ActsMatrix<eBoundSize, eAlignmentSize>
-
using AlignmentToPathMatrix = ActsMatrix<1, eAlignmentSize>
-
using AlignmentToPositionMatrix = ActsMatrix<3, eAlignmentSize>
-
using AlignmentVector = ActsVector<eAlignmentSize>
-
using AngleAxis3 = Eigen::AngleAxis<ActsScalar>
-
using Any = AnyBase<sizeof(void*)>
-
using AxisScalar = Vector3::Scalar
-
using BoundaryIntersection = ObjectIntersection<BoundarySurface, Surface>
Intersection with a
BoundarySurface
.
-
using BoundarySurface = BoundarySurfaceT<TrackingVolume>
BoundarySurface of a volume.
-
using BoundarySurfacePtr = std::shared_ptr<const BoundarySurfaceT<AbstractVolume>>
-
using BoundMatrix = ActsMatrix<eBoundSize, eBoundSize>
-
using BoundSquareMatrix = ActsSquareMatrix<eBoundSize>
-
using BoundToFreeMatrix = ActsMatrix<eFreeSize, eBoundSize>
-
using BoundTrackParameters = GenericBoundTrackParameters<SinglyCharged>
-
using BoundVariantMeasurement = VariantMeasurement<BoundIndices>
Variant that can hold all possible bound measurements.
-
using BoundVector = ActsVector<eBoundSize>
-
using CalibrationContext = std::any
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
-
using ColorRGB = std::array<int, 3>
The color typedef.
It’s an array of three numbers [0, 255] representing the RGB color values.
-
template<typename T>
using ConstTrackAccessor = TrackAccessorBase<T, true>
-
using CurvilinearTrackParameters = GenericCurvilinearTrackParameters<SinglyCharged>
-
using DefaultExtension = detail::GenericDefaultExtension<double>
A typedef for the default GenericDefaultExtension with double.
-
using DenseEnvironmentExtension = detail::GenericDenseEnvironmentExtension<double>
A typedef for the default GenericDenseEnvironmentExtension with double.
-
using EAxis = Acts::detail::EquidistantAxis
-
using Envelope = std::array<ActsScalar, 2>
-
using FreeMatrix = ActsMatrix<eFreeSize, eFreeSize>
-
using FreeSquareMatrix = ActsSquareMatrix<eFreeSize>
-
using FreeToBoundMatrix = ActsMatrix<eBoundSize, eFreeSize>
-
using FreeToPathMatrix = ActsMatrix<1, eFreeSize>
-
using FreeTrackParameters = GenericFreeTrackParameters<SinglyCharged>
-
using FreeVariantMeasurement = VariantMeasurement<FreeIndices>
Variant that can hold all possible free measurements.
-
using FreeVector = ActsVector<eFreeSize>
-
using GeometryContext = ContextType
This is the central definition of the Acts payload object regarding detector geometry status (e.g.
alignment)
It is propagated through the code to allow for event/thread dependent geometry changes
-
using Grid2D = Acts::detail::Grid<Acts::AccumulatedVolumeMaterial, EAxis, EAxis>
-
using HashedString = std::uint32_t
-
using Intersection2D = Intersection<2>
-
using Intersection3D = Intersection<3>
-
typedef BinnedArray<LayerPtr> LayerArray
A BinnedArray to a std::shared_ptr of a layer.
Layers are constructed with shared_ptr factories, hence the layer array is describes as:
-
using LayerIntersection = ObjectIntersection<Layer, Surface>
-
typedef std::shared_ptr<const Layer> LayerPtr
A std::shared_ptr to a Layer.
-
typedef std::vector<LayerPtr> LayerVector
A vector of std::shared_ptr to layers.
-
using MagneticFieldContext = ContextType
This is the central definition of the Acts payload object regarding magnetic field status.
It is propagated through the code to allow for event/thread dependent magnetic field changes
-
using MaterialGridAxisData = std::tuple<double, double, size_t>
-
using MaterialSlabMatrix = std::vector<MaterialSlabVector>
-
using MaterialSlabVector = std::vector<MaterialSlab>
-
typedef std::shared_ptr<Layer> MutableLayerPtr
-
typedef std::shared_ptr<TrackingVolume> MutableTrackingVolumePtr
-
typedef std::vector<MutableTrackingVolumePtr> MutableTrackingVolumeVector
-
using NetworkBatchInput = Eigen::Array<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
-
using NeutralBoundTrackParameters = GenericBoundTrackParameters<Neutral>
-
using NeutralCurvilinearTrackParameters = GenericCurvilinearTrackParameters<Neutral>
-
using NeutralFreeTrackParameters = GenericFreeTrackParameters<Neutral>
-
using OrientedSurfaces = std::vector<OrientedSurface>
-
using ProjectorBitset = uint64_t
-
using Ray3D = Ray<ActsScalar, 3>
-
using RecordedMaterialTrack = std::pair<std::pair<Acts::Vector3, Acts::Vector3>, RecordedMaterial>
And recorded material track.
this is start: position, start momentum and the Recorded material
-
using RecordedMaterialVolumePoint = std::vector<std::pair<Acts::MaterialSlab, std::vector<Acts::Vector3>>>
list of point used in the mapping of a volume
-
using RotationMatrix2 = ActsMatrix<2, 2>
-
using RotationMatrix3 = ActsMatrix<3, 3>
-
template<typename source_link_iterator_t>
using SourceLinkAccessorDelegate = Delegate<std::pair<source_link_iterator_t, source_link_iterator_t>(const Surface&)> Delegate type that retrieves a range of source links to for a given surface to be processed by the CKF.
-
using SourceLinkSurfaceAccessor = Delegate<const Surface*(const SourceLink&)>
Delegate to unpack the surface associated with a source link.
-
template<typename external_spacepoint_t>
using SpacePointGrid = detail::Grid<std::vector<std::unique_ptr<InternalSpacePoint<external_spacepoint_t>>>, detail::Axis<detail::AxisType::Equidistant, detail::AxisBoundaryType::Closed>, detail::Axis<detail::AxisType::Variable, detail::AxisBoundaryType::Bound>>
-
using SquareMatrix2 = ActsSquareMatrix<2>
-
using SquareMatrix3 = ActsSquareMatrix<3>
-
using SquareMatrix4 = ActsSquareMatrix<4>
-
typedef std::tuple<std::shared_ptr<const Acts::Surface>, std::shared_ptr<const Acts::ISurfaceMaterial>, Acts::GeometryContext> SurfaceAndMaterialWithContext
-
typedef ObjectIntersection<Surface> SurfaceIntersection
Typedef of the surface intersection.
-
using SurfaceMatcher = std::function<bool(const GeometryContext &gctx, BinningValue, const Surface*, const Surface*)>
-
using surfaceMaterialPointer = const Acts::ISurfaceMaterial*
-
using SurfaceMatrix = std::vector<SurfaceVector>
-
typedef std::shared_ptr<const Surface> SurfacePtr
-
typedef std::vector<SurfacePtr> SurfacePtrVector
-
typedef std::vector<const Surface*> SurfaceVector
-
template<typename T>
using TrackAccessor = TrackAccessorBase<T, false>
-
using TrackIndexType = std::uint32_t
-
typedef std::pair<const Acts::TrackingVolume*, std::shared_ptr<const Acts::IVolumeMaterial>> TrackingVolumeAndMaterial
-
typedef BinnedArray<TrackingVolumePtr> TrackingVolumeArray
A BinnedArray of a std::shared_tr to a TrackingVolume.
-
using TrackingVolumeBoundaries = std::vector<TrackingVolumeBoundaryPtr>
-
using TrackingVolumeBoundaryPtr = std::shared_ptr<const BoundarySurfaceT<TrackingVolume>>
-
using TrackingVolumeOrderPosition = std::pair<TrackingVolumePtr, Vector3>
-
using TrackingVolumePointer = const Acts::TrackingVolume*
-
typedef std::shared_ptr<const TrackingVolume> TrackingVolumePtr
A std::shared_ptr to a tracking volume.
-
typedef std::vector<TrackingVolumePtr> TrackingVolumeVector
A std::vector of a std::shared_ptr to a TrackingVolume.
-
using Transform2 = Eigen::Transform<ActsScalar, 2, Eigen::AffineCompact>
-
using Transform3 = Eigen::Transform<ActsScalar, 3, Eigen::Affine>
-
using Translation2 = Eigen::Translation<ActsScalar, 2>
-
using Translation3 = Eigen::Translation<ActsScalar, 3>
-
using V3Matrix = std::vector<V3Vector>
-
using V3Vector = std::vector<Vector3>
-
using VariantCovariance = std::variant<BoundSquareMatrix, FreeSquareMatrix>
-
template<typename indices_t>
using VariantMeasurement = typename detail::VariantMeasurementGenerator<indices_t, detail::kParametersSize<indices_t>>::Type Variant that can contain all possible measurements in a parameter space.
- Template Parameters
indices_t – Parameter index type, determines the full parameter space
-
using VariantTransportJacobian = std::variant<BoundMatrix, BoundToFreeMatrix, FreeToBoundMatrix, FreeMatrix>
-
using Vector2 = ActsVector<2>
-
using Vector3 = ActsVector<3>
-
using Vector4 = ActsVector<4>
-
typedef std::shared_ptr<const VolumeBounds> VolumeBoundsPtr
-
using volumeMaterialPointer = const Acts::IVolumeMaterial*
Enums
-
enum AlignmentIndices
Components of alignment parameters vector.
To be used to access components by named indices instead of just numbers. This must be a regular
enum
and not a scopedenum class
to allow implicit conversion to an integer. The enum value are thus visible directly innamespace Acts
and are prefixed to avoid naming collisions.Values:
-
enumerator eAlignmentCenter0
-
enumerator eAlignmentCenter1
-
enumerator eAlignmentCenter2
-
enumerator eAlignmentRotation0
-
enumerator eAlignmentRotation1
-
enumerator eAlignmentRotation2
-
enumerator eAlignmentSize
-
enumerator eAlignmentCenter0
-
enum BinningOption
flag for open/closed bins
Values:
-
enumerator open
-
enumerator closed
-
enumerator open
-
enum BinningType
, BinningOption & BinningAccess
BinningType:
Enumeration to qualify the binning type for the use of the LayerArrayCreator and the TrackingVolumeArrayCreator
BinningOption: open: [0,max] closed: 0 -> nextbin -> max -> 0
BinningValue necessary access to global positions
Values:
-
enumerator equidistant
-
enumerator arbitrary
-
enum BinningValue
how to take the global / local position
Values:
-
enumerator binX
-
enumerator binY
-
enumerator binZ
-
enumerator binR
-
enumerator binPhi
-
enumerator binRPhi
-
enumerator binH
-
enumerator binEta
-
enumerator binMag
-
enumerator binValues
-
enumerator binX
-
enum BoundarySurfaceFace
Enum to describe the position of the BoundarySurface respectively to the frame orientatin of the volume, this is mainly meant for code readability.
The different numeration sequences can be found in the documentation of the actual VolumeBounds implementations.
The order of faces is chosen to follow - as much as possible - a circular structure.
Values:
-
enumerator negativeFaceXY
-
enumerator positiveFaceXY
-
enumerator negativeFaceYZ
-
enumerator positiveFaceYZ
-
enumerator negativeFaceZX
-
enumerator positiveFaceZX
-
enumerator cylinderCover
-
enumerator tubeInnerCover
-
enumerator tubeOuterCover
-
enumerator tubeSectorNegativePhi
-
enumerator tubeSectorPositivePhi
-
enumerator tubeSectorInnerCover
-
enumerator tubeSectorOuterCover
-
enumerator trapezoidFaceAlpha
-
enumerator trapezoidFaceBeta
-
enumerator index0
-
enumerator index1
-
enumerator index2
-
enumerator index3
-
enumerator index4
-
enumerator index5
-
enumerator index6
-
enumerator index7
-
enumerator index8
-
enumerator index9
-
enumerator index10
-
enumerator index11
-
enumerator undefinedFace
-
enumerator negativeFaceXY
-
enum BoundIndices
Components of a bound track parameters vector.
To be used to access components by named indices instead of just numbers. This must be a regular
enum
and not a scopedenum class
to allow implicit conversion to an integer. The enum value are thus visible directly innamespace Acts
and are prefixed to avoid naming collisions.Values:
-
enumerator eBoundLoc0
-
enumerator eBoundLoc1
-
enumerator eBoundPhi
-
enumerator eBoundTheta
-
enumerator eBoundQOverP
-
enumerator eBoundTime
-
enumerator eBoundSize
-
enumerator eBoundLoc0
-
enum class CombinatorialKalmanFilterError
Values:
-
enumerator UpdateFailed
-
enumerator SmoothFailed
-
enumerator OutputConversionFailed
-
enumerator MeasurementSelectionFailed
-
enumerator PropagationReachesMaxSteps
-
enumerator UpdateFailed
-
enum CoordinateIndices
Components of coordinate vectors.
To be used to access coordinate components by named indices instead of magic numbers. This must be a regular
enum
and not a scopedenum class
to allow implicit conversion to an integer. The enum value are thus visible directly innamespace Acts
.This index enum is not user-configurable (in contrast e.g. to the track parameter index enums) since it must be compatible with varying dimensionality (2d-4d) and other access methods (
.{x,y,z}()
accessors).Values:
-
enumerator ePos0
-
enumerator ePos1
-
enumerator ePos2
-
enumerator eTime
-
enumerator eMom0
-
enumerator eMom1
-
enumerator eMom2
-
enumerator eEnergy
-
enumerator eX
-
enumerator eY
-
enumerator eZ
-
enumerator ePos0
-
enum class DelegateType
Ownership enum for
Delegate
.Values:
-
enumerator Owning
-
enumerator NonOwning
-
enumerator Owning
-
enum class DetectorMeasurementInfo : short
Values:
-
enumerator eDefault
-
enumerator eDetailed
-
enumerator eDefault
-
enum class EigenStepperError
Values:
-
enumerator StepSizeStalled
-
enumerator StepInvalid
-
enumerator StepSizeAdjustmentFailed
-
enumerator StepSizeStalled
-
enum FreeIndices
Components of a free track parameters vector.
To be used to access components by named indices instead of just numbers. This must be a regular
enum
and not a scopedenum class
to allow implicit conversion to an integer. The enum value are thus visible directly innamespace Acts
and are prefixed to avoid naming collisions.Values:
-
enumerator eFreePos0
-
enumerator eFreePos1
-
enumerator eFreePos2
-
enumerator eFreeTime
-
enumerator eFreeDir0
-
enumerator eFreeDir1
-
enumerator eFreeDir2
-
enumerator eFreeQOverP
-
enumerator eFreeSize
-
enumerator eFreePos0
-
enum class IntersectionStatus : int
Status enum.
Values:
-
enumerator missed
-
enumerator unreachable
-
enumerator reachable
-
enumerator onSurface
-
enumerator missed
-
enum class KalmanFitterError
Values:
-
enumerator ForwardUpdateFailed
-
enumerator BackwardUpdateFailed
-
enumerator SmoothFailed
-
enumerator OutputConversionFailed
-
enumerator NoMeasurementFound
-
enumerator ForwardUpdateFailed
-
enum LayerType
For code readability, it distinguishes between different type of layers, which steers the behaviour in the navigation.
Values:
-
enumerator passive
-
enumerator active
-
enumerator passive
-
enum LinIndices
Enum to access the components of a track parameter vector.
Here, we parametrize the track via a 4D point on the track, the momentum angles of the particle at that point, and q/p or 1/p.
Note
It would make sense to rename these parameters if they are used outside of track linearization.
Note
This must be a regular
enum
and not a scopedenum class
to allow implicit conversion to an integer. The enum value are thus visible directly innamespace Acts
and are prefixed to avoid naming collisions.Values:
-
enumerator eLinPos0
-
enumerator eLinPos1
-
enumerator eLinPos2
-
enumerator eLinTime
-
enumerator eLinPhi
-
enumerator eLinTheta
-
enumerator eLinQOverP
-
enumerator eLinSize
-
enumerator eLinPosSize
-
enumerator eLinMomSize
-
enumerator eLinPos0
-
enum class MagneticFieldError
Values:
-
enumerator OutOfBounds
-
enumerator OutOfBounds
-
enum MappingType
This enum describes the type of surface material mapping.
Values:
-
enumerator PreMapping
-
enumerator Default
-
enumerator PostMapping
-
enumerator Sensor
-
enumerator PreMapping
-
enum class MaterialUpdateStage : int
This is a steering enum to tell which material update stage:
PreUpdate : update on approach of a surface
FullUpdate : update when passing a surface
PostUpdate : update when leaving a surface
Values:
-
enumerator PreUpdate
-
enumerator FullUpdate
-
enumerator PostUpdate
-
enum class MixtureReductionMethod
Available reduction methods for the reduction of a Gaussian mixture.
Values:
-
enumerator eMean
-
enumerator eMaxWeight
-
enumerator eMean
-
enum class MultiStepperError
Values:
-
enumerator ComponentNotOnSurface
-
enumerator StateOfMultipleComponentsRequested
-
enumerator AverageTrackLeftCurrentVolume
-
enumerator AllComponentsSteppingError
-
enumerator AllComponentsConversionToBoundFailed
-
enumerator SomeComponentsConversionToBoundFailed
-
enumerator ComponentNotOnSurface
-
enum NoiseUpdateMode
to tell how to deal with noise term in covariance transport
removeNoise: subtract noise term
addNoise: add noise term
Values:
-
enumerator removeNoise
-
enumerator addNoise
-
enum PdgParticle
Symbolic values for commonly used PDG particle numbers.
Values:
-
enumerator eInvalid
-
enumerator eElectron
-
enumerator eAntiElectron
-
enumerator ePositron
-
enumerator eMuon
-
enumerator eAntiMuon
-
enumerator eTau
-
enumerator eAntiTau
-
enumerator eGamma
-
enumerator ePionZero
-
enumerator ePionPlus
-
enumerator ePionMinus
-
enumerator eNeutron
-
enumerator eAntiNeutron
-
enumerator eProton
-
enumerator eAntiProton
-
enumerator eLead
-
enumerator eInvalid
-
enum class PropagatorError
Values:
-
enumerator Failure
-
enumerator WrongDirection
-
enumerator StepCountLimitReached
-
enumerator Failure
-
enum class SpacePointCandidateType : short
Values:
-
enumerator eBottom
-
enumerator eTop
-
enumerator eBottom
-
enum class SurfaceError
Values:
-
enumerator GlobalPositionNotOnSurface
-
enumerator GlobalPositionNotOnSurface
-
enum TrackStateFlag
This enum describes the type of TrackState.
Values:
-
enumerator MeasurementFlag
-
enumerator ParameterFlag
-
enumerator OutlierFlag
-
enumerator HoleFlag
-
enumerator MaterialFlag
-
enumerator NumTrackStateFlags
-
enumerator MeasurementFlag
-
enum class TrackStatePropMask : std::uint8_t
Collection of bit masks to enable steering which components of a track state should be initialized, and which should be left invalid.
These mask values can be combined using binary operators, so (TrackStatePropMask::Predicted | TrackStatePropMask::Jacobian) will instruct allocating storage for both predicted parameters (including covariance) and a jacobian. The enum is used as a strong type wrapper around the bits to prevent autoconversion from integer
Values:
-
enumerator None
-
enumerator Predicted
-
enumerator Filtered
-
enumerator Smoothed
-
enumerator Jacobian
-
enumerator Calibrated
-
enumerator All
-
enumerator None
-
enum class VertexingError
Values:
-
enumerator NumericFailure
-
enumerator EmptyInput
-
enumerator SeedingError
-
enumerator NotConverged
-
enumerator ElementNotFound
-
enumerator NoCovariance
-
enumerator NumericFailure
-
enum WrappingCondition
Values:
-
enumerator Undefined
inconsistency detected
-
enumerator Attaching
attach the volumes
-
enumerator Inserting
insert the new volume
-
enumerator Wrapping
wrap the new volume around
-
enumerator CentralInserting
insert the new one into the center
-
enumerator CentralWrapping
wrap the new central volume around
-
enumerator NoWrapping
no inner volume present - no wrapping needed
-
enumerator Undefined
Functions
-
ACTS_STATIC_CHECK_CONCEPT(ChargeConcept, AnyCharge)
-
ACTS_STATIC_CHECK_CONCEPT(ChargeConcept, Neutral)
-
ACTS_STATIC_CHECK_CONCEPT(ChargeConcept, NonNeutralCharge)
-
ACTS_STATIC_CHECK_CONCEPT(ChargeConcept, SinglyCharged)
-
ACTS_STATIC_CHECK_CONCEPT(ConstMultiTrajectoryBackend, ConstVectorMultiTrajectory)
-
ACTS_STATIC_CHECK_CONCEPT(ConstTrackContainerBackend, ConstVectorTrackContainer)
-
ACTS_STATIC_CHECK_CONCEPT(MutableMultiTrajectoryBackend, VectorMultiTrajectory)
-
ACTS_STATIC_CHECK_CONCEPT(TrackContainerBackend, VectorTrackContainer)
-
void addCylinderLayerProtoMaterial(dd4hep::DetElement detElement, Layer &cylinderLayer, const Logger &logger = getDummyLogger())
Helper method to translate DD4hep material to Acts::ISurfaceMaterial.
This is used to assign proto material to Cylinder Layers
- Parameters
detElement – the DD4hep detector element for which this material is assigned
cylinderLayer – is the target layer
logger – a
Logger
for output
-
void addDiscLayerProtoMaterial(dd4hep::DetElement detElement, Layer &discLayer, const Logger &logger = getDummyLogger())
Helper method to translate DD4hep material to Acts::ISurfaceMaterial.
Thisis used to assign proto material to Disc Layers
- Parameters
detElement – the DD4hep detector element for which this material is assigned
discLayer – is the target layer
logger – a
Logger
for output
-
void addLayerProtoMaterial(const dd4hep::rec::VariantParameters ¶ms, Layer &layer, const std::vector<std::pair<const std::string, Acts::BinningOption>> &binning, const Logger &logger = getDummyLogger())
Helper method to be called for Cylinder and Disc Proto material.
For both, cylinder and disc, the closed binning value is “binPhi”
- Parameters
params – An instance of
DD4hep::VariantParameters
layer – the Layer to assign the proto material
binning – the Binning prescription for the ActsExtension
logger – a
Logger
for output
-
BinUtility adjustBinUtility(const BinUtility &bu, const CuboidVolumeBounds &cBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of cuboid volume bounds
- Parameters
bu – BinUtility at source
cBounds – the Cuboid volume bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const CutoutCylinderVolumeBounds &cBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of cutout cylinder volume bounds
- Parameters
bu – BinUtility at source
cBounds – the Cutout Cylinder volume bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const CylinderBounds &cBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of cylinder bounds
- Parameters
bu – BinUtility at source
cBounds – the Cylinder bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const CylinderVolumeBounds &cBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of cylinder volume bounds
- Parameters
bu – BinUtility at source
cBounds – the Cylinder volume bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const RadialBounds &rBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of radial bounds
- Parameters
bu – BinUtility at source
rBounds – the Radial bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const RectangleBounds &pBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of plane bounds
- Parameters
bu – BinUtility at source
pBounds – the Rectangle bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const Surface &surface, const GeometryContext &gctx)
adjust the BinUtility bu to a surface
- Parameters
bu – BinUtility at source
surface – Surface to which the adjustment is being done
gctx – Geometry context to get the surfaces transform
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const TrapezoidBounds &pBounds, const Transform3 &transform)
adjust the BinUtility bu to the dimensions of plane bounds
- Parameters
bu – BinUtility at source
pBounds – the Trapezoid bounds to adjust to
transform – Transform for the adjusted
BinUtility
- Returns
new updated BinUtiltiy
-
BinUtility adjustBinUtility(const BinUtility &bu, const Volume &volume)
adjust the BinUtility bu to a volume
- Parameters
bu – BinUtility at source
volume – Volume to which the adjustment is being done
- Returns
new updated BinUtiltiy
-
inline const std::vector<std::string> &binningValueNames()
screen output option
-
template<typename MatrixType>
MatrixType bitsetToMatrix(const std::bitset<MatrixType::RowsAtCompileTime * MatrixType::ColsAtCompileTime> bs) Convert a bitset to a matrix of integers, with each element set to the bit value.
Note
How the bits are assigned to matrix elements depends on the storage type of the matrix being converted (row-major or col-major)
- Template Parameters
MatrixType – Matrix type that is produced
- Parameters
bs – The bitset to convert
- Returns
A matrix with the integer values of the bits from
bs
-
template<typename A, typename B>
inline ActsMatrix<A::RowsAtCompileTime, B::ColsAtCompileTime> blockedMult(const A &a, const B &b) Perform a blocked matrix multiplication, avoiding Eigen GEMM methods.
- Template Parameters
A – The type of the first matrix, which should be MxN
B – The type of the second matrix, which should be NxP
- Parameters
a – [in] An MxN matrix of type A
b – [in] An NxP matrix of type P
- Returns
The product ab
-
template<typename track_container_t, typename track_state_container_t, template<typename> class holder_t>
void calculateTrackQuantities(Acts::TrackProxy<track_container_t, track_state_container_t, holder_t, false> track) Helper function to calculate a number of track level quantities and store them on the track itself.
Note
The input track needs to be mutable, so
ReadOnly=false
- Template Parameters
track_container_t – the track container backend
track_state_container_t – the track state container backend
holder_t – the holder type for the track container backends
- Parameters
track – A mutable track proxy to operate on
-
template<typename T, typename U>
T clampValue(U value) Clamp a numeric value to another type, respecting range of the target type.
- Template Parameters
T – the target type
U – the source type
- Parameters
value – the value to clamp
- Returns
the clamped value
-
void collectCompounds_dd4hep(dd4hep::DetElement &detElement, std::vector<dd4hep::DetElement> &compounds)
Method internally used by convertDD4hepDetector to collect all volumes of a compound detector.
- Parameters
detElement – [in] the dd4hep::DetElement of the volume of which the compounds should be collected
compounds – [out] the DD4hep::DetElements of the compounds contained by detElement
-
void collectLayers_dd4hep(dd4hep::DetElement &detElement, std::vector<dd4hep::DetElement> &layers, const Logger &logger)
Method internally used by convertDD4hepDetector.
- Parameters
detElement – [in] the dd4hep::DetElement of the volume of which the layers should be collected
layers – [out] the DD4hep::DetElements of the layers contained by detElement
logger – a
Logger
for output
-
void collectSubDetectors_dd4hep(dd4hep::DetElement &detElement, std::vector<dd4hep::DetElement> &subdetectors, const Logger &logger)
Method internally used by convertDD4hepDetector to collect all sub detectors Sub detector means each ‘compound’ DetElement or DetElements which are declared as ‘isBarrel’ or ‘isBeampipe’ by their extension.
- Parameters
detElement – [in] the dd4hep::DetElement of the volume of which the sub detectors should be collected
subdetectors – [out] the DD4hep::DetElements of the sub detectors contained by detElement
logger – a
Logger
for output
-
float computeEnergyLossBethe(const MaterialSlab &slab, float m, float qOverP, float absQ)
Compute the mean energy loss due to ionisation and excitation.
This computes the mean energy loss -dE(x) through a material with the given properties, i.e. it computes
where -dE/dx is given by the Bethe formula. The computations are valid for intermediate particle energies.-dE(x) = -dE/dx * x
- Parameters
slab – The traversed material and its properties
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float computeEnergyLossLandau(const MaterialSlab &slab, float m, float qOverP, float absQ)
Compute the most propable energy loss due to ionisation and excitation.
Compute the mean energy loss due to ionisation and excitation.
This computes the mean energy loss -dE(x) through a material with the given properties, i.e. it computes
where -dE/dx is given by the Bethe formula. The computations are valid for intermediate particle energies.-dE(x) = -dE/dx * x
This computes the most probable energy loss -dE(x) through a material of the given properties and thickness as described by the mode of the Landau-Vavilov-Bichsel distribution. The computations are valid for intermediate particle energies.
- Parameters
slab – The traversed material and its properties
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float computeEnergyLossLandauFwhm(const MaterialSlab &slab, float m, float qOverP, float absQ)
Compute the full with half maximum of landau energy loss distribution.
See also
computeEnergyLossBethe for parameters description
-
float computeEnergyLossLandauSigma(const MaterialSlab &slab, float m, float qOverP, float absQ)
Compute the Gaussian-equivalent sigma for the ionisation loss fluctuations.
This is the sigma parameter of a Gaussian distribution with the same full-width-half-maximum as the Landau-Vavilov-Bichsel distribution. The computations are valid for intermediate particle energies.
See also
computeEnergyLossBethe for parameters description
-
float computeEnergyLossLandauSigmaQOverP(const MaterialSlab &slab, float m, float qOverP, float absQ)
Compute q/p Gaussian-equivalent sigma due to ionisation loss fluctuations.
Compute the mean energy loss due to ionisation and excitation.
This computes the mean energy loss -dE(x) through a material with the given properties, i.e. it computes
where -dE/dx is given by the Bethe formula. The computations are valid for intermediate particle energies.-dE(x) = -dE/dx * x
- Parameters
slab – The traversed material and its properties
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float computeEnergyLossMean(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Compute the combined mean energy loss.
This computes the combined mean energy loss -dE(x) including ionisation and radiative effects. The computations are valid over a wide range of particle energies.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float computeEnergyLossMode(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Compute the combined most probably energy loss.
Compute the combined mean energy loss.
This computes the combined mean energy loss -dE(x) including ionisation and radiative effects. The computations are valid over a wide range of particle energies.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float computeEnergyLossRadiative(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Compute the mean energy loss due to radiative effects at high energies.
This computes the mean energy loss -dE(x) using an approximative formula. Bremsstrahlung is always included; direct e+e- pair production and photo-nuclear interactions only for muons.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float computeMultipleScatteringTheta0(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Compute the core width of the projected planar scattering distribution.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
Global method which creates the TrackingGeometry from DD4hep input.
This method returns a std::unique_ptr of the TrackingGeometry from the World DD4hep DetElement.
- Attention
The default thickness should be set thin enough that no touching or overlapping with the next layer can happen.
Note
Possible binningtypes:
arbitrary - of the sizes if the surfaces and the distance in between vary. This mode finds out the bin boundaries by scanning through the surfaces.
equidistant - if the sensitive surfaces are placed equidistantly
Note
equidistant binningtype is recommended because it is faster not only while building the geometry but also for look up during the extrapolation
Note
Layers containing surfaces per default are not allowed to be attached to each other (navigation will fail at this point). However, to allow material layers (not containing surfaces) to be attached to each other, this default thickness is needed. In this way, the layer will be thin (with space to the next layer), but the material will have the’real’ thickness.
- Parameters
worldDetElement – [in] the DD4hep DetElement of the world
logger – [in] A logger instance geometry building
bTypePhi – [in] is how the sensitive surfaces (modules) should be binned in a layer in phi direction.
bTypeR – [in] is how the sensitive surfaces (modules) should be binned in a layer in r direction
bTypeZ – [in] is how the sensitive surfaces (modules) should be binned in a layer in z direction
layerEnvelopeR – [in] the tolerance added to the geometrical extension in r of the layers contained to build the volume envelope around
layerEnvelopeZ – [in] the tolerance added to the geometrical extension in z of the layers contained to build the volume envelope around
defaultLayerThickness – [in] In case no surfaces (to be contained by the layer) are handed over, the layer thickness will be set to this value
sortSubDetectors – [in]
std::function
which should be used in order to sort all sub detectors (=all Detelements collected by the methodcollectSubDetectors()
) from bottom to top to ensure correct wrapping of the volumes, which is needed for navigation. Therefore the different hierarchies need to be sorted ascending. The default is sorting by ID.gctx – The geometry context to use
matDecorator – is the material decorator that loads material maps
geometryIdentifierHook – Hook to apply to surfaces during geometry closure.
- Throws
std::logic_error – if an error in the translation occurs
- Pre
Before using this method make sure, that the preconditions described in DD4hepPlugins are met.
- Returns
std::unique_ptr to the full TrackingGeometry
The Tracking geometry needs to be built from bottom to top to ensure Navigation. Therefore the different hierarchies need to be sorted ascending. Per default the sub detectors are sorted by the id of their dd4hep::DetElement. In case another sorting needs to be applied, the users can provide their own function
-
Grid2D createGrid(MaterialGridAxisData gridAxis1, MaterialGridAxisData gridAxis2)
Helper method that creates the cache grid for the mapping.
This grid allows the collection of material at a the anchor points.
Note
The data of the axes is given in the std::array as {minimum value, maximum value, number of bins}
- Parameters
gridAxis1 – [in] Axis data
gridAxis2 – [in] Axis data
- Returns
The grid
-
Grid3D createGrid(MaterialGridAxisData gridAxis1, MaterialGridAxisData gridAxis2, MaterialGridAxisData gridAxis3)
Helper method that creates the cache grid for the mapping.
This grid allows the collection of material at a the anchor points.
Note
The data of the axes is given in the std::array as {minimum value, maximum value, number of bins}
- Parameters
gridAxis1 – [in] Axis data
gridAxis2 – [in] Axis data
gridAxis3 – [in] Axis data
- Returns
The grid
-
Grid2D createGrid2D(const BinUtility &bins, std::function<Acts::Vector2(Acts::Vector3)> &transfoGlobalToLocal)
Create a 2DGrid using a BinUtility.
Also determine the corresponding global to local transform and grid mapping function
- Parameters
bins – [in] BinUtility of the volume to be mapped
transfoGlobalToLocal – [in] Global to local transform to be updated.
- Returns
the 3D grid
-
Grid3D createGrid3D(const BinUtility &bins, std::function<Acts::Vector3(Acts::Vector3)> &transfoGlobalToLocal)
Create a 3DGrid using a BinUtility.
Also determine the corresponding global to local transform and grid mapping function
- Parameters
bins – [in] BinUtility of the volume to be mapped
transfoGlobalToLocal – [in] Global to local transform to be updated.
- Returns
the 3D grid
-
std::shared_ptr<Acts::ProtoSurfaceMaterial> createProtoMaterial(const dd4hep::rec::VariantParameters ¶ms, const std::string &valueTag, const std::vector<std::pair<const std::string, Acts::BinningOption>> &binning, const Logger &logger = getDummyLogger())
Helper method to create proto material - to be called from the addProto(…) methods.
- Parameters
params – An instance of
DD4hep::VariantParameters
valueTag – the xml tag for to ActsExtension to be parsed
binning – the Binning prescription for the ActsExtension
logger – a
Logger
for output
-
std::shared_ptr<const CylinderVolumeHelper> cylinderVolumeHelper_dd4hep(const Logger &logger)
Helper method internally used to create a default Acts::CylinderVolumeBuilder.
-
template<class T, class decorator_t>
void decorateJson([[maybe_unused]] const decorator_t *decorator, [[maybe_unused]] const T &src, [[maybe_unused]] nlohmann::json &dest)
-
template<class T, class decorator_t>
void decorateJson([[maybe_unused]] const decorator_t *decorator, [[maybe_unused]] const T *src, [[maybe_unused]] nlohmann::json &dest)
-
float deriveEnergyLossBetheQOverP(const MaterialSlab &slab, float m, float qOverP, float absQ)
Derivative of the Bethe energy loss with respect to q/p.
Compute the mean energy loss due to ionisation and excitation.
This computes the mean energy loss -dE(x) through a material with the given properties, i.e. it computes
where -dE/dx is given by the Bethe formula. The computations are valid for intermediate particle energies.-dE(x) = -dE/dx * x
- Parameters
slab – The traversed material and its properties
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float deriveEnergyLossLandauQOverP(const MaterialSlab &slab, float m, float qOverP, float absQ)
Derivative of the most probable ionisation energy loss with respect to q/p.
Compute the mean energy loss due to ionisation and excitation.
This computes the mean energy loss -dE(x) through a material with the given properties, i.e. it computes
where -dE/dx is given by the Bethe formula. The computations are valid for intermediate particle energies.-dE(x) = -dE/dx * x
- Parameters
slab – The traversed material and its properties
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float deriveEnergyLossMeanQOverP(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Derivative of the combined mean energy loss with respect to q/p.
Compute the combined mean energy loss.
This computes the combined mean energy loss -dE(x) including ionisation and radiative effects. The computations are valid over a wide range of particle energies.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float deriveEnergyLossModeQOverP(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Derivative of the combined most probable energy loss with respect to q/p.
Compute the combined mean energy loss.
This computes the combined mean energy loss -dE(x) including ionisation and radiative effects. The computations are valid over a wide range of particle energies.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
float deriveEnergyLossRadiativeQOverP(const MaterialSlab &slab, PdgParticle absPdg, float m, float qOverP, float absQ)
Derivative of the mean radiative energy loss with respect to q/p.
Compute the mean energy loss due to radiative effects at high energies.
This computes the mean energy loss -dE(x) using an approximative formula. Bremsstrahlung is always included; direct e+e- pair production and photo-nuclear interactions only for muons.
- Parameters
slab – The traversed material and its properties
absPdg – Absolute particle type PDG identifier
m – Particle mass
qOverP – Particle charge divided by absolute momentum
absQ – Absolute particle charge
-
template<typename container_type, typename container_type_iter = decltype(std::begin(std::declval<container_type>())), typename = decltype(std::end(std::declval<container_type>()))>
constexpr auto enumerate(container_type &&iterable) Helper utility to allow indexed enumeration with structured binding.
Usage:
for (auto [ i, value ] = enumerate(container) ) { … };
with ‘container’ any stl-like container
-
template<typename spacepoint_iterator_t>
std::optional<BoundVector> estimateTrackParamsFromSeed(const GeometryContext &gctx, spacepoint_iterator_t spBegin, spacepoint_iterator_t spEnd, const Surface &surface, const Vector3 &bField, ActsScalar bFieldMin, const Acts::Logger &logger = getDummyLogger(), ActsScalar mass = 139.57018 * UnitConstants::MeV) Estimate the full track parameters from three space points.
This method is based on the conformal map transformation. It estimates the full bound track parameters, i.e. (loc0, loc1, phi, theta, q/p, t) at the bottom space point. The bottom space is assumed to be the first element in the range defined by the iterators. It must lie on the surface provided for the representation of the bound track parameters. The magnetic field (which might be along any direction) is also necessary for the momentum estimation.
It resembles the method used in ATLAS for the track parameters estimated from seed, i.e. the function InDet::SiTrackMaker_xk::getAtaPlane here: https://acode-browser.usatlas.bnl.gov/lxr/source/athena/InnerDetector/InDetRecTools/SiTrackMakerTool_xk/src/SiTrackMaker_xk.cxx
- Template Parameters
spacepoint_iterator_t – The type of space point iterator
- Parameters
gctx – is the geometry context
spBegin – is the begin iterator for the space points
spEnd – is the end iterator for the space points
surface – is the surface of the bottom space point. The estimated bound track parameters will be represented also at this surface
bField – is the magnetic field vector
bFieldMin – is the minimum magnetic field required to trigger the estimation of q/pt
logger – A logger instance
mass – is the estimated particle mass
- Returns
optional bound parameters
-
template<typename spacepoint_iterator_t>
std::optional<BoundVector> estimateTrackParamsFromSeed(spacepoint_iterator_t spBegin, spacepoint_iterator_t spEnd, const Logger &logger = getDummyLogger()) Estimate the track parameters on the xy plane from at least three space points.
It assumes the trajectory projection on the xy plane is a circle, i.e. the magnetic field is along global z-axis.
The method is based on V. Karimaki NIM A305 (1991) 187-191: https://doi.org/10.1016/0168-9002(91)90533-V
no weights are used in Karimaki’s fit; d0 is the distance of the point of closest approach to the origin, 1/R is the curvature, phi is the angle of the direction propagation (counter clockwise as positive) at the point of cloest approach.
- Template Parameters
spacepoint_iterator_t – The type of space point iterator
- Parameters
spBegin – is the begin iterator for the space points
spEnd – is the end iterator for the space points
logger – A logger instance
- Returns
optional bound track parameters with the estimated d0, phi and 1/R stored with the indices, eBoundLoc0, eBoundPhi and eBoundQOverP, respectively. The bound parameters with other indices are set to zero.
-
Acts::InterpolatedBFieldMap<Acts::detail::Grid<Acts::Vector2, Acts::detail::EquidistantAxis, Acts::detail::EquidistantAxis>> fieldMapRZ(const std::function<size_t(std::array<size_t, 2> binsRZ, std::array<size_t, 2> nBinsRZ)> &localToGlobalBin, std::vector<double> rPos, std::vector<double> zPos, std::vector<Acts::Vector2> bField, double lengthUnit = UnitConstants::mm, double BFieldUnit = UnitConstants::T, bool firstQuadrant = false)
Method to setup the FieldMap.
e.g.: we have small grid with the values: r={2,3}, z ={4,5}, the corresponding indices are i (belonging to r) and j (belonging to z), the globalIndex is M (belonging to the value of the magnetic field B(r,z)) and the field map is:
r
i
z
j
B(r,z)
M
2
0
4
0
2.323
0
2
0
5
1
2.334
1
3
1
4
0
2.325
2
3
1
5
1
2.331
3
In this case the function would look like:
[](std::array<size_t, 2> binsRZ, std::array<size_t, 2> nBinsRZ) { return (binsRZ.at(0) * nBinsRZ.at(1) + binsRZ.at(1)); }
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The function localToGlobalBin determines how the magnetic field was stored in the vector in respect to the grid values
Note
If
firstQuadrant
is true will create a field that is symmetric for all quadrants. e.g. we have the grid values r={0,1} with BFieldValues={2,3} on the r axis. If the flag is set to true the r-axis grid values will be set to {-1,0,1} and the BFieldValues will be set to {3,2,3}.- Parameters
localToGlobalBin – Function mapping the local bins of r,z to the global bin of the map magnetic field value
rPos – [in] Values of the grid points in r
zPos – [in] Values of the grid points in z
bField – [in] The magnetic field values inr r and z for all given grid points stored in a vector
lengthUnit – [in] The unit of the grid points
BFieldUnit – [in] The unit of the magnetic field
firstQuadrant – [in] Flag if set to true indicating that only the first quadrant of the grid points and the BField values has been given.
- Returns
A field map instance for use in interpolation.
-
Acts::InterpolatedBFieldMap<Acts::detail::Grid<Acts::Vector3, Acts::detail::EquidistantAxis, Acts::detail::EquidistantAxis, Acts::detail::EquidistantAxis>> fieldMapXYZ(const std::function<size_t(std::array<size_t, 3> binsXYZ, std::array<size_t, 3> nBinsXYZ)> &localToGlobalBin, std::vector<double> xPos, std::vector<double> yPos, std::vector<double> zPos, std::vector<Acts::Vector3> bField, double lengthUnit = UnitConstants::mm, double BFieldUnit = UnitConstants::T, bool firstOctant = false)
Method to setup the FieldMap.
e.g.: we have small grid with the values: x={2,3}, y={3,4}, z ={4,5}, the corresponding indices are i (belonging to x), j (belonging to y) and k (belonging to z), the globalIndex is M (belonging to the value of the magnetic field B(x,y,z)) and the field map is:
x
i
y
j
z
k
B(x,y,z)
M
2
0
3
0
4
0
2.323
0
2
0
3
0
5
1
2.334
1
2
0
4
1
4
0
2.325
2
2
0
4
1
5
1
2.331
3
3
1
3
0
4
0
2.323
4
3
1
3
0
5
1
2.334
5
3
1
4
1
4
0
2.325
6
3
1
4
1
5
1
2.331
7
In this case the function would look like:
[](std::array<size_t, 3> binsXYZ, std::array<size_t, 3> nBinsXYZ) { return (binsXYZ.at(0) * (nBinsXYZ.at(1) * nBinsXYZ.at(2)) + binsXYZ.at(1) * nBinsXYZ.at(2) + binsXYZ.at(2)); }
Note
The grid point values
xPos
,yPos
andzPos
do not need to be sorted or unique (this will be done inside the function)Note
The function
localToGlobalBin
determines how the magnetic field was stored in the vector in respect to the grid valuesNote
If
firstOctant
is true, the function will assume a symmetrical field for all quadrants. e.g. we have the grid values z={0,1} with BFieldValues={2,3} on the r axis. If the flag is set to true the z-axis grid values will be set to {-1,0,1} and the BFieldValues will be set to {3,2,3}.- Parameters
localToGlobalBin – Function mapping the local bins of x,y,z to the global bin of the map magnetic field value
xPos – [in] Values of the grid points in x
yPos – [in] Values of the grid points in y
zPos – [in] Values of the grid points in z
bField – [in] The magnetic field values inr r and z for all given grid points stored in a vector
lengthUnit – [in] The unit of the grid points
BFieldUnit – [in] The unit of the magnetic field
firstOctant – [in] Flag if set to true indicating that only the first octant of the grid points and the BField values has been given.
- Returns
A field map instance for use in interpolation.
-
std::optional<float> findCharge(PdgParticle pdg)
Find the charge for a given PDG particle number.
- Returns
Charge in native units.
-
std::optional<float> findMass(PdgParticle pdg)
Find the mass for a given PDG particle number.
- Returns
Mass in native units.
-
std::optional<std::string_view> findName(PdgParticle pdg)
Find a descriptive particle name for a given PDG particle number.
- Returns
Particle name.
-
std::optional<ParticleData> findParticleData(PdgParticle pdg)
Find all known particle data for a given PDG particle number.
- Returns
Particle name.
-
void from_json(const nlohmann::json &j, BinningData &bd)
-
void from_json(const nlohmann::json &j, BinUtility &bu)
-
void from_json(const nlohmann::json &j, Extent &e)
-
void from_json(const nlohmann::json &j, Material &t)
-
void from_json(const nlohmann::json &j, MaterialSlab &t)
-
void from_json(const nlohmann::json &j, MaterialSlabMatrix &t)
-
void from_json(const nlohmann::json &j, ProtoDetector &pd)
-
void from_json(const nlohmann::json &j, ProtoVolume &pv)
-
void from_json(const nlohmann::json &j, surfaceMaterialPointer &material)
-
void from_json(const nlohmann::json &j, Transform3 &t)
-
void from_json(const nlohmann::json &j, volumeMaterialPointer &material)
-
std::unique_ptr<const Logger> getDefaultLogger(const std::string &name, const Logging::Level &lvl, std::ostream *log_stream = &std::cout)
get default debug output logger
This function returns a pointer to a Logger instance with the following decorations enabled:
time stamps
name of logging instance
debug level
- Parameters
name – [in] name of the logger instance
lvl – [in] debug threshold level
log_stream – [in] output stream used for printing debug messages
- Returns
pointer to logging instance
-
const Logger &getDummyLogger()
-
template<typename T>
T getParam(const std::string &key, dd4hep::DetElement &elt) Helper function to extract a parameter value from a dd4hep detector element from VariantParameters.
- Template Parameters
T – The value type
- Parameters
key – The key of the value to extract
elt – The detector element instance
- Returns
A copy of the value contained in the params instance
-
template<typename T>
T getParamOr(const std::string &key, const dd4hep::DetElement &elt, T alternative) Get a parameter value or an alternative value if either the VariantParameters extension isn’t set, or it doesn’t contain the demanded key.
- Template Parameters
T – The value type
- Parameters
key – The key of the value to extract
elt – The detector element instance
alternative – The value to return if no params are set of the key doesn’t exist
- Returns
The value behind key, or
alternative
-
inline const dd4hep::rec::VariantParameters &getParams(const dd4hep::DetElement &elt)
Helper function to extract a VariantParameters instance, const version.
- Parameters
elt – The detector element instance
- Returns
The VariantParameters instance
-
inline dd4hep::rec::VariantParameters &getParams(dd4hep::DetElement &elt)
Helper function to extract a VariantParameters instance.
- Parameters
elt – The detector element instance
- Returns
The VariantParameters instance
-
std::function<double(Acts::Vector3)> globalToLocalFromBin(Acts::BinningValue &type)
return a function that return the coordinate corresponding to type of bin
- Parameters
type – [in] Type of bin
- Returns
a coordinate transform function
-
constexpr HashedString hashString(std::string_view s)
-
inline bool hasParam(const std::string &key, dd4hep::DetElement &elt)
Check if a detector element has a key set in its VariantParameters.
- Parameters
key – The key to check existence for
elt – The detector element instance
- Returns
True if the element has VariantParameters and the key exists, false if either of these is not true
-
inline bool hasParams(dd4hep::DetElement &elt)
Check if a detector element has VariantParameters set.
- Parameters
elt – The detector element instance
- Returns
True if the VariantParameters exist, false if not
-
template<typename T, size_t N, class Point1, class Point2 = Point1, class Point3 = Point2, typename = std::enable_if_t<detail::can_interpolate<Point1, Point2, Point3, T>::value>>
inline T interpolate(const Point1 &position, const Point2 &lowerCorner, const Point3 &upperCorner, const std::array<T, N> &values) performs linear interpolation inside a hyper box
Note
Given
U
andV
of value typeT
as well as twodouble
a
andb
, then the following must be a valid expressiona * U + b * V
yielding an object which is (implicitly) convertible toT
.All
Point
types must represent d-dimensional positions and support coordinate access usingoperator
[] which should return adouble
(or a value which is implicitly convertible). Coordinate indices must start at 0.N
is the number of hyper box corners which is \(2^d\) where \(d\) is the dimensionality of the hyper box. The dimensionality must be consistent with the providedPoint
types.Definition of the canonical order for sorting the field values: The hyper box corners are numbered according to the following scheme. Each corner is defined by the set of lower/upper boundary limits in each dimension
i
. This can be represented by a binary code (from left to right) where a0
stands for a lower bound along this axis and a1
stand for the upper bound along this axis. The left most bit corresponds to the first dimension and the bits to the left correspond to the 2nd, 3rd… dimension. The binary code can be interpreted as integer which gives the number of the corresponding hyper box corner. The field values are ordered according to ascending hyper box corner numbers.
As an example assume we have a 3D box with
lowerCorner
= (1,2,3) andupperCorner
= (4,5,6). The eight corners with their bit patterns and corner numbers are:(1,2,3): 000 = 0
(1,2,6): 001 = 1
(1,5,3): 010 = 2
(1,5,6): 011 = 3
(4,2,3): 100 = 4
(4,2,6): 101 = 5
(4,5,3): 110 = 6
(4,5,6): 111 = 7
- Template Parameters
T – type of values to be interpolated
N – number of hyper box corners
Point1 – type specifying geometric positions
Point2 – type specifying geometric positions
Point3 – type specifying geometric positions
- Parameters
position – [in] position to which to interpolate
lowerCorner – [in] generalized lower-left corner of hyper box (containing the minima of the hyper box along each dimension)
upperCorner – [in] generalized upper-right corner of hyper box (containing the maxima of the hyper box along each dimension)
values – [in] field values at the hyper box corners sorted in the canonical order defined below.
- Returns
interpolated value at given position
- Pre
position
must describe a position inside the given hyper box, that is \(\text{lowerCorner}[i] \le \text{position}[i] \le \text{upperCorner}[i] \quad \forall i=0, \dots, d-1\).
-
std::error_code make_error_code(Acts::CombinatorialKalmanFilterError e)
-
std::error_code make_error_code(Acts::KalmanFitterError e)
-
std::error_code make_error_code(Acts::MagneticFieldError e)
-
std::error_code make_error_code(Acts::MultiStepperError e)
-
std::error_code make_error_code(Acts::PropagatorError e)
-
std::error_code make_error_code(Acts::SurfaceError e)
-
std::error_code make_error_code(Acts::VertexingError e)
-
std::error_code make_error_code(EigenStepperError e)
-
template<typename box_t>
box_t *make_octree(std::vector<std::unique_ptr<box_t>> &store, const std::vector<box_t*> &prims, size_t max_depth = 1, typename box_t::value_type envelope1 = 0) Build an octree from a list of bounding boxes.
Note
store
andprims
do not need to contain the same objects.store
is only used to pass ownership back to the caller while preserving memory location.- Template Parameters
box_t – Works will all box types.
- Parameters
store – Owns the created boxes by means of
std::unique_ptr
.prims – Boxes to store. This is a read only vector.
max_depth – No subdivisions beyond this level.
envelope1 – Envelope to add/subtract to dimensions in all directions.
- Returns
Pointer to the top most bounding box, containing the entire octree
-
static inline constexpr PdgParticle makeAbsolutePdgParticle(PdgParticle pdg)
Convert an anti-particle to its particle and leave particles as-is.
-
template<typename InputVector>
inline auto makeCurvilinearUnitU(const Eigen::MatrixBase<InputVector> &direction) Construct the first curvilinear unit vector
U
for the given direction.The special case of the direction vector pointing along the z-axis is handled by forcing the unit vector to along the x-axis.
- Parameters
direction – is the input direction vector
- Returns
a normalized vector in the x-y plane orthogonal to the direction.
-
template<typename InputVector>
inline auto makeCurvilinearUnitVectors(const Eigen::MatrixBase<InputVector> &direction) Construct the curvilinear unit vectors
U
andV
for the given direction.With
T
the normalized input direction, the three vectorsU
,V
, andT
form an orthonormal basis set, i.e. they satisfywith the additional condition thatU x V = T V x T = U T x U = V
U
is located in the global x-y plane.- Parameters
direction – is the input direction vector
- Returns
normalized unit vectors
U
andV
orthogonal to the direction.
-
template<typename T>
inline Eigen::Matrix<T, 3, 1> makeDirectionFromPhiEta(T phi, T eta) Construct a normalized direction vector from phi angle and pseudorapidity.
Note
The input arguments intentionally use the same template type so that a compile error occurs if inconsistent input types are used. Avoids unexpected implicit type conversions and forces the user to explicitly cast mismatched input types.
- Parameters
phi – is the direction angle in the x-y plane.
eta – is the pseudorapidity towards the z-axis.
-
template<typename T>
inline Eigen::Matrix<T, 3, 1> makeDirectionFromPhiTheta(T phi, T theta) Construct a normalized direction vector from phi and theta angle.
Note
The input arguments intentionally use the same template type so that a compile error occurs if inconsistent input types are used. Avoids unexpected implicit type conversions and forces the user to explicitly cast mismatched input types.
- Parameters
phi – is the direction angle in radian in the x-y plane.
theta – is the polar angle in radian towards the z-axis.
-
template<typename parameters_t, typename covariance_t, typename indices_t, typename ...tail_indices_t>
auto makeMeasurement(SourceLink source, const Eigen::MatrixBase<parameters_t> ¶ms, const Eigen::MatrixBase<covariance_t> &cov, indices_t index0, tail_indices_t... tailIndices) -> Measurement<indices_t, 1u + sizeof...(tail_indices_t)> Construct a fixed-size measurement for the given indices.
This helper function can be used to create a fixed-size measurement using an explicit set of indices, e.g.
for a 2d measurement w/ one position and time.auto m = makeMeasurement(s, p, c, eBoundLoc0, eBoundTime);
Note
The indices must be ordered and must be consistent with the content of parameters and covariance.
- Template Parameters
parameters_t – Input parameters vector type
covariance_t – Input covariance matrix type
indices_t – Parameter index type, determines the full parameter space
tail_indices_t – Helper types required to support variadic arguments; all types must be convertibale to
indices_t
.
- Parameters
source – The link that connects to the underlying detector readout
params – Measured parameters values
cov – Measured parameters covariance
index0 – Required parameter index, measurement must be at least 1d
tailIndices – Additional parameter indices for larger measurements
- Returns
Fixed-size measurement w/ the correct type and given inputs
-
template<typename T>
inline Eigen::Matrix<T, 2, 1> makePhiThetaFromDirection(Eigen::Matrix<T, 3, 1> unitDir) Construct a phi and theta angle from a direction vector.
- Parameters
unitDir – 3D vector indicating a direction
-
MaterialGrid2D mapMaterialPoints(Grid2D &grid)
Average the material collected in a 2D grid and use it to create a 2D material grid.
- Parameters
grid – [in] The material collecting grid coordinate
- Returns
The average material grid decomposed into classification numbers
-
MaterialGrid3D mapMaterialPoints(Grid3D &grid)
Average the material collected in a 3D grid and use it to create a 3D material grid.
- Parameters
grid – [in] The material collecting grid coordinate
- Returns
The average material grid decomposed into classification numbers
-
MaterialMapper<detail::Grid<Material::ParametersVector, detail::EquidistantAxis, detail::EquidistantAxis>> materialMapperRZ(const std::function<size_t(std::array<size_t, 2> binsRZ, std::array<size_t, 2> nBinsRZ)> &materialVectorToGridMapper, std::vector<double> rPos, std::vector<double> zPos, const std::vector<Acts::Material> &material, double lengthUnit = UnitConstants::mm)
Method to setup the MaterialMapper.
e.g.: we have small grid with the values: r={2,3}, z ={4,5}, the corresponding indices are i (belonging to r) and j (belonging to z), the globalIndex is M (belonging to the values of the Material) and the map is:
r
i
z
j
M
2
0
4
0
0
2
0
5
1
1
3
1
4
0
2
3
1
5
1
3
In this case the function would look like:
[](std::array<size_t, 2> binsRZ, std::array<size_t, 2> nBinsRZ) { return (binsRZ.at(0) * nBinsRZ.at(1) + binsRZ.at(1)); }
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The function localToGlobalBin determines how the material was stored in the vector in respect to the grid values
- Parameters
materialVectorToGridMapper – [in] Function mapping the vector of material to the map of material values
rPos – [in] Values of the grid points in r
zPos – [in] Values of the grid points in z
material – [in] The material classification values in r and z for all given grid points stored in a vector
lengthUnit – [in] The unit of the grid points
-
MaterialMapper<detail::Grid<Material::ParametersVector, detail::EquidistantAxis, detail::EquidistantAxis, detail::EquidistantAxis>> materialMapperXYZ(const std::function<size_t(std::array<size_t, 3> binsXYZ, std::array<size_t, 3> nBinsXYZ)> &materialVectorToGridMapper, std::vector<double> xPos, std::vector<double> yPos, std::vector<double> zPos, const std::vector<Material> &material, double lengthUnit = UnitConstants::mm)
Method to setup the MaterialMapper.
e.g.: we have small grid with the values: x={2,3}, y={3,4}, z ={4,5}, the corresponding indices are i (belonging to x), j (belonging to y) and k (belonging to z), the globalIndex is M (belonging to the values of the Material) and the map is:
x
i
y
j
z
k
M
2
0
3
0
4
0
0
2
0
3
0
5
1
1
2
0
4
1
4
0
2
2
0
4
1
5
1
3
3
1
3
0
4
0
4
3
1
3
0
5
1
5
3
1
4
1
4
0
6
3
1
4
1
5
1
7
In this case the function would look like:
[](std::array<size_t, 3> binsXYZ, std::array<size_t, 3> nBinsXYZ) { return (binsXYZ.at(0) * (nBinsXYZ.at(1) * nBinsXYZ.at(2)) + binsXYZ.at(1) * nBinsXYZ.at(2) + binsXYZ.at(2)); }
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The values do not need to be sorted or unique (this will be done inside the function)
Note
The function localToGlobalBin determines how the material was stored in the vector in respect to the grid values
- Parameters
materialVectorToGridMapper – [in] Function mapping the vector of material to the map of material values
xPos – [in] Values of the grid points in x
yPos – [in] Values of the grid points in y
zPos – [in] Values of the grid points in z
material – [in] The material classification values in x, y and z for all given grid points stored in a vector
lengthUnit – [in] The unit of the grid points
-
template<typename Derived>
auto matrixToBitset(const Eigen::PlainObjectBase<Derived> &m) Convert an integer matrix to a bitset.
Note
How the bits are ordered depends on the storage type of the matrix being converted (row-major or col-major)
- Template Parameters
Derived – Eigen base concrete type
- Parameters
m – Matrix that is converted
- Returns
The converted bitset.
-
template<typename T>
std::array<typename T::value_type, 2u> min_max(const T &tseries) Return min/max from a (optionally) sorted series, obsolete with C++20 (ranges)
- Template Parameters
T – a numeric series
- Parameters
tseries – is the number series
- Returns
[ min, max ] in an array of length 2
-
inline bool operator!=(const DigitizationSourceLink &lhs, const DigitizationSourceLink &rhs)
-
inline bool operator!=(const SurfaceBounds &lhs, const SurfaceBounds &rhs)
-
inline constexpr double operator*(Direction dir, double value)
-
inline constexpr float operator*(Direction dir, float value)
-
inline constexpr int operator*(Direction dir, int value)
-
inline constexpr double operator*(double value, Direction dir)
-
inline constexpr float operator*(float value, Direction dir)
-
inline constexpr int operator*(int value, Direction dir)
-
inline constexpr double operator*=(double &value, Direction dir)
-
inline constexpr float operator*=(float &value, Direction dir)
-
inline constexpr int operator*=(int &value, Direction dir)
-
inline std::ostream &operator<<(std::ostream &os, BoundarySurfaceFace &face)
-
template<typename T, typename U, size_t V>
std::ostream &operator<<(std::ostream &os, const AxisAlignedBoundingBox<T, U, V> &box) Overload of the << operator for bounding boxes.
- Template Parameters
T – entity type
U – value type
V – dimension
- Parameters
os – The output stream
box – The bounding box
- Returns
The given output stream.
-
inline std::ostream &operator<<(std::ostream &os, const ConstrainedStep &step)
-
inline std::ostream &operator<<(std::ostream &os, const IVisualization3D &hlp)
Overload of the << operator to facilitate writing to streams.
- Parameters
os – The output stream
hlp – The helper instance
-
std::ostream &operator<<(std::ostream &os, const Material &material)
-
std::ostream &operator<<(std::ostream &os, const MaterialSlab &materialSlab)
-
template<typename T, size_t D>
std::ostream &operator<<(std::ostream &os, const Ray<T, D> &ray) Overload of the outstream operator.
- Parameters
os – The out stream
ray – The ray to write to
os
- Returns
The outstream given in
os
-
inline std::ostream &operator<<(std::ostream &os, const std::tuple<const Surface&, const GeometryContext&> &tup)
Print surface information to the provided stream.
Internally invokes the
surface.toStream(...)
-method. This can be easily used e.g. likestd::cout << std::tie(surface, geometryContext);
-
inline std::ostream &operator<<(std::ostream &os, const SurfaceBounds &sb)
-
inline std::ostream &operator<<(std::ostream &os, const TrackSelector::Config &cuts)
- Parameters
os – Output stream
cuts – Cuts to print
- Returns
Reference to the output stream
-
inline std::ostream &operator<<(std::ostream &os, const TrackSelector::EtaBinnedConfig &cfg)
- Parameters
os – Output stream
cfg – Configuration to print
- Returns
Reference to the output stream
-
template<typename indices_t>
std::ostream &operator<<(std::ostream &os, const VariantMeasurement<indices_t> &vm)
-
std::ostream &operator<<(std::ostream &os, Direction dir)
-
std::ostream &operator<<(std::ostream &os, GeometryIdentifier id)
-
inline std::ostream &operator<<(std::ostream &os, IntersectionStatus status)
Ostream-operator for the IntersectionStatus enum.
-
std::ostream &operator<<(std::ostream &os, MaterialUpdateStage matUpdate)
-
std::ostream &operator<<(std::ostream &os, PdgParticle pdg)
Print PDG particle numbers with a descriptive name.
-
std::ostream &operator<<(std::ostream &os, TrackStatePropMask mask)
-
std::ostream &operator<<(std::ostream &sl, const BinUtility &bgen)
Overload of << operator for std::ostream for debug output.
-
std::ostream &operator<<(std::ostream &sl, const Extent &rhs)
Overload of << operator for std::ostream for debug output.
-
std::ostream &operator<<(std::ostream &sl, const GlueVolumesDescriptor &gvd)
-
std::ostream &operator<<(std::ostream &sl, const Volume &vol)
Overload of << operator for std::ostream for debug output.
-
std::ostream &operator<<(std::ostream &sl, const VolumeBounds &vb)
Overload of << operator for std::ostream for debug output.
-
inline bool operator==(const DigitizationSourceLink &lhs, const DigitizationSourceLink &rhs)
-
inline bool operator==(const SurfaceBounds &lhs, const SurfaceBounds &rhs)
-
inline bool operator==(const VolumeBounds &lhs, const VolumeBounds &rhs)
-
std::optional<std::string_view> pdgToShortAbsString(PdgParticle pdg)
-
template<typename T>
std::tuple<typename T::value_type, ActsScalar> range_medium(const T &tseries) Return range and medium of a sorted numeric series.
- Template Parameters
T – a numeric series
- Parameters
tseries – is the number series
- Returns
[ range, medium ] in an tuple
-
template<typename mixture_t, typename projector_t = Acts::Identity>
auto reduceGaussianMixture(const mixture_t &mixture, const Surface &surface, MixtureReductionMethod method, projector_t &&projector = projector_t{}) Reduce Gaussian mixture.
- Parameters
mixture – The mixture iterable
surface – The surface, on which the bound state is
method – How to reduce the mixture
projector – How to project a element of the iterable to something like a std::tuple< double, BoundVector, BoundMatrix >
- Returns
parameters and covariance as std::tuple< BoundVector, BoundMatrix >
-
template<typename T>
constexpr T safeExp(T val) noexcept Calculate the exponential function while avoiding FPEs.
- Parameters
val – argument for which the exponential function should be evaluated.
- Returns
0 in the case of underflow, std::numeric_limits<T>::infinity in the case of overflow, std::exp(val) else
-
template<typename MatrixType, typename ResultType = MatrixType>
std::optional<ResultType> safeInverse(const MatrixType &m) noexcept FPE “safe” functions.
Our main motivation for this is that users might have a strict FPE policy which would flag every single occurrence as a failure and then sombody has to investigate. Since we are processing a high number of events and floating point numbers sometimes work in mysterious ways the caller of this function might want to hide FPEs and handle them in a more controlled way. Calculate the inverse of an Eigen matrix after checking if it can be numerically inverted. This allows to catch potential FPEs before they occur.
- Template Parameters
Derived – Eigen derived concrete type
Result – Eigen result type defaulted to input type
- Parameters
m – Eigen matrix to invert
- Returns
The theta value
-
Acts::InterpolatedBFieldMap<Acts::detail::Grid<Acts::Vector2, Acts::detail::EquidistantAxis, Acts::detail::EquidistantAxis>> solenoidFieldMap(std::pair<double, double> rlim, std::pair<double, double> zlim, std::pair<size_t, size_t> nbins, const SolenoidBField &field)
Function which takes an existing SolenoidBField instance and creates a field mapper by sampling grid points from the analytical solenoid field.
- Parameters
rlim – pair of r bounds
zlim – pair of z bounds
nbins – pair of bin counts
field – the solenoid field instance
- Returns
A field map instance for use in interpolation.
-
inline void sortDetElementsByID(std::vector<dd4hep::DetElement> &det)
Sort function which sorts dd4hep::DetElement by their ID.
- Parameters
det – [inout] the dd4hep::DetElements to be sorted
-
std::shared_ptr<Surface> surfaceFromJson(const nlohmann::json &j)
Conversion to Surface from jsonn.
- Parameters
j – the read-in json object
- Returns
a shared_ptr to a surface object for type polymorphism
Conversion to Surface from json in correct type.
The type is given as a template argument in order to be able to construct the correct fitting types for surfaces.
- Parameters
j – the read-in json object
- Returns
a shared_ptr to a typed surface object for type polymorphism
-
template<template<size_t> class Callable, size_t N, size_t NMAX, typename ...Args>
auto template_switch(size_t v, Args&&... args) Dispatch a call based on a runtime value on a function taking the value at compile time.
This function allows to write a templated functor, which accepts a
size_t
like parameter at compile time. It is then possible to make a call to the corresponding instance of the functor based on a runtime value. To achieve this, the function essentially created a if cascade betweenN
andNMAX
, attempting to find the right instance. Because the cascade is visible to the compiler entirely, it should be able to optimize.Note
Callable
is expected to have a static member functioninvoke
that is callable withArgs
- Template Parameters
Callable – Type which takes a size_t as a compile time param
N – Value from which to start the dispatch chain, i.e. 0 in most cases
NMAX – Maximum value up to which to attempt a dispatch
- Parameters
v – The runtime value to dispatch on
args – Additional arguments passed to
Callable::invoke()
.
-
template<size_t N, size_t NMAX, typename Lambda, typename ...Args>
auto template_switch_lambda(size_t v, Lambda &&func, Args&&... args) Alternative version of
template_switch
which accepts a generic lambda and communicates the dimension via an integral constant type.- Template Parameters
N – Value from which to start the dispatch chain, i.e. 0 in most cases
NMAX – Maximum value up to which to attempt a dispatch
- Parameters
v – The runtime value to dispatch on
func – The lambda to invoke
args – Additional arguments passed to
func
-
template<std::size_t kDIM, typename value_type>
std::array<value_type, kDIM> to_array(const std::vector<value_type> &vecvals) This can be abandoned with C++20 to use the std::to_array method.
Note
only the first kDIM elements will obviously be filled, if the vector tends to be longer, it is truncated
- Parameters
vecvals – the vector of bound values to be converted
- Returns
an array with the filled values
-
void to_json(nlohmann::json &j, const BinningData &bd)
-
void to_json(nlohmann::json &j, const BinUtility &bu)
-
void to_json(nlohmann::json &j, const Extent &e)
-
void to_json(nlohmann::json &j, const Material &t)
-
void to_json(nlohmann::json &j, const MaterialSlab &t)
-
void to_json(nlohmann::json &j, const ProtoDetector &pd)
-
void to_json(nlohmann::json &j, const ProtoVolume &pv)
Non-contextual conversion of a surface.
Note
it will take the default context
-
void to_json(nlohmann::json &j, const Surface &surface)
Non-contextual conversion of a surface.
Note
it will take the default context
-
void to_json(nlohmann::json &j, const SurfaceAndMaterialWithContext &surface)
Conversion of a pair of surface and material used for the material mapping.
-
void to_json(nlohmann::json &j, const SurfaceBounds &bounds)
-
void to_json(nlohmann::json &j, const surfaceMaterialPointer &material)
-
void to_json(nlohmann::json &j, const TrackingVolumeAndMaterial &volume)
Conversion of a pair of tracking volume and material used for the material mapping.
-
void to_json(nlohmann::json &j, const TrackingVolumePointer &volume)
Conversion of a const pointer on a tracking volume used to write the geometry.
Conversion of a tracking volume.
-
void to_json(nlohmann::json &j, const Transform3 &t)
-
void to_json(nlohmann::json &j, const VolumeBounds &bounds)
-
void to_json(nlohmann::json &j, const volumeMaterialPointer &material)
Contextual conversion of a surface.
- Parameters
j – the json to be filled
surface – the surface to be converted
gctx – the geometry context for this
-
inline std::string toString(const Acts::Transform3 &transform, int precision = 4, const std::string &offset = "")
Print out a transform in a structured way.
- Parameters
transform – The transform to print
precision – Numeric output precision
offset – Offset in front of matrix lines
- Returns
The printed string
-
inline std::string toString(const Acts::Translation3 &translation, int precision = 4)
Print out a translation in a structured way.
- Parameters
translation – The translation to print
precision – Numeric output precision
- Returns
The printed string
-
template<typename derived_t>
inline std::string toString(const Eigen::MatrixBase<derived_t> &matrix, int precision = 4, const std::string &offset = "") Print out a matrix in a structured way.
- Template Parameters
derived_t – Type of the matrix
- Parameters
matrix – The matrix to print
precision – Numeric output precision
offset – Offset in front of matrix lines
- Returns
The printed string
-
template<ACTS_CONCEPT track_container_t(TrackContainerBackend), typename traj_t>
TrackContainer(const track_container_t &container, const traj_t &traj) -> TrackContainer<track_container_t, traj_t, detail::ConstRefHolder>
-
template<ACTS_CONCEPT track_container_t(TrackContainerBackend), typename traj_t>
TrackContainer(track_container_t &&container, traj_t &&traj) -> TrackContainer<track_container_t, traj_t, detail::ValueHolder>
-
template<ACTS_CONCEPT track_container_t(TrackContainerBackend), typename traj_t>
TrackContainer(track_container_t &container, traj_t &traj) -> TrackContainer<track_container_t, traj_t, detail::RefHolder>
-
template<typename external_spacepoint_t, typename callable_t>
void transformCoordinates(Acts::SpacePointData &spacePointData, const std::vector<external_spacepoint_t*> &vec, const external_spacepoint_t &spM, bool bottom, std::vector<LinCircle> &linCircleVec, callable_t &&extractFunction)
-
template<typename external_spacepoint_t>
void transformCoordinates(Acts::SpacePointData &spacePointData, const std::vector<InternalSpacePoint<external_spacepoint_t>*> &vec, const InternalSpacePoint<external_spacepoint_t> &spM, bool bottom, std::vector<LinCircle> &linCircleVec) Transform a vector of spacepoints to u-v space circles with respect to a given middle spacepoint.
- Template Parameters
external_spacepoint_t – The external spacepoint type.
- Parameters
spacePointData – [in] Auxiliary variables used by the seeding
vec – [in] The list of bottom or top spacepoints
spM – [in] The middle spacepoint.
bottom – [in] Should be true if vec are bottom spacepoints.
linCircleVec – [out] The output vector to write to.
-
template<typename external_spacepoint_t, typename callable_t>
LinCircle transformCoordinates(const external_spacepoint_t &sp, const external_spacepoint_t &spM, bool bottom, callable_t &&extractFunction)
-
template<typename external_spacepoint_t>
LinCircle transformCoordinates(const InternalSpacePoint<external_spacepoint_t> &sp, const InternalSpacePoint<external_spacepoint_t> &spM, bool bottom) Transform two spacepoints to a u-v space circle.
This function is a non-vectorized version of transformCoordinates.
- Template Parameters
external_spacepoint_t – The external spacepoint type.
- Parameters
sp – [in] The first spacepoint to use, either a bottom or top.
spM – [in] The middle spacepoint to use.
bottom – [in] Should be true if sp is a bottom SP.
-
std::tuple<VariantCovariance, VariantTransportJacobian> transportCovarianceToBound(const GeometryContext &gctx, const Surface &surface, const FreeVector ¶meters, CovarianceCache &cCache)
Transport the covariance to a new (surface) bound definition.
- Parameters
gctx – [in] The current geometry context
surface – [in] The surface of the bound representation
parameters – [in] The free parametrisation on the surface
cCache – [in,out] the covariance cache (to be modified)
- Returns
a variant transport jacobian
-
std::tuple<VariantCovariance, VariantTransportJacobian> transportCovarianceToCurvilinear(const Vector3 &direction, CovarianceCache &cCache)
Transport the covariance to a new (curvilinear) bound definition.
- Parameters
direction – [in] The track direction at the new curvilinear definition
cCache – [in,out] the covariance cache (to be modified)
-
std::tuple<VariantCovariance, VariantTransportJacobian> transportCovarianceToFree(CovarianceCache &cCache)
Transport the covariance to a new free definition.
- Parameters
cCache – [in,out] the covariance cache (to be modified)
Helper function to unpack a vector of
shared_ptr
into a vector of raw pointers.- Template Parameters
T – the stored type
- Parameters
items – The vector of
shared_ptr
- Returns
The unpacked vector
Helper function to unpack a vector of
shared_ptr
into a vector of raw pointers (const version)- Template Parameters
T – the stored type
- Parameters
items – The vector of
shared_ptr
- Returns
The unpacked vector
Helper function to unpack a vector of
shared_ptr
into a vector of raw pointers.- Template Parameters
T – the stored type
- Parameters
items – The vector of
shared_ptr
- Returns
The unpacked vector
-
template<typename L, typename A, typename B>
auto visit_measurement(A &¶m, B &&cov, size_t dim, L &&lambda) Dispatch a lambda call on an overallocated parameter vector and covariance matrix, based on a runtime dimension value.
Inside the lambda call, the vector and matrix will have fixed dimensions, but will still point back to the originally given overallocated values.
Note
No requirements on
A
andB
are made, to enable a single overload for both const and non-const matrices/vectors.- Template Parameters
L – The lambda type
A – The parameter vector type
B – The covariance matrix type
- Parameters
param – The parameter vector
cov – The covariance matrix
dim – The actual dimension as a runtime value
lambda – The lambda to call with the statically sized subsets
-
template<typename L>
auto visit_measurement(size_t dim, L &&lambda) Dispatch a generic lambda on a measurement dimension.
This overload doesn’t assume anything about what is needed inside the lambda, it communicates the dimension via an integral constant type
- Template Parameters
L – The generic lambda type to call
- Parameters
dim – The runtime dimension of the measurement
lambda – The generic lambda instance to call
- Returns
Returns the lambda return value
-
std::shared_ptr<const CylinderVolumeBuilder> volumeBuilder_dd4hep(dd4hep::DetElement subDetector, const Logger &logger, BinningType bTypePhi = equidistant, BinningType bTypeR = equidistant, BinningType bTypeZ = equidistant, double layerEnvelopeR = UnitConstants::mm, double layerEnvelopeZ = UnitConstants::mm, double defaultLayerThickness = UnitConstants::fm)
Method internally used to create an Acts::CylinderVolumeBuilder.
This method creates an Acts::CylinderVolumeBuilder from a sub detector (= ‘compound’ DetElement or DetElements which are declared as ‘isBarrel’ or ‘isBeampipe’ by their extension.
- Attention
The default thickness should be set thin enough that no touching or overlapping with the next layer can happen.
Note
Possible binningtypes:
arbitrary - of the sizes if the surfaces and the distance in between vary. This mode finds out the bin boundaries by scanning through the surfaces.
equidistant - if the sensitive surfaces are placed equidistantly
Note
equidistant binningtype is recommended because it is faster not only while building the geometry but also for look up during the extrapolation
Note
Layers containing surfaces per default are not allowed to be attached to each other (navigation will fail at this point). However, to allow material layers (not containing surfaces) to be attached to each other, this default thickness is needed. In this way, the layer will be thin (with space to the next layer), but the material will have the’real’ thickness.
- Parameters
subDetector – [in] the DD4hep DetElement of the subdetector
logger – [in] A logger instance
bTypePhi – [in] is how the sensitive surfaces (modules) should be binned in a layer in phi direction.
bTypeR – [in] is how the sensitive surfaces (modules) should be binned in a layer in r direction
bTypeZ – [in] is how the sensitive surfaces (modules) should be binned in a layer in z direction
layerEnvelopeR – [in] the tolerance added to the geometrical extension in r of the layers contained to build the volume envelope around
layerEnvelopeZ – [in] the tolerance added to the geometrical extension in z of the layers contained to build the volume envelope around
defaultLayerThickness – [in] In case no surfaces (to be contained by the layer) are handed over, the layer thickness will be set to this value
- Returns
std::shared_ptr the Acts::CylinderVolumeBuilder which can be used to build the full tracking geometry
-
template<typename external_spacepoint_t>
bool xyzCoordinateCheck(Acts::SpacePointData &spacePointData, const Acts::SeedFinderConfig<external_spacepoint_t> &config, const Acts::InternalSpacePoint<external_spacepoint_t> &sp, const double *spacepointPosition, double *outputCoordinates) Check the compatibility of spacepoint coordinates in xyz assuming the Bottom-Middle direction with the strip meassument details.
- Template Parameters
external_spacepoint_t – The external spacepoint type.
- Parameters
spacePointData – [in] Auxiliary variables used by the seeding
config – [in] SeedFinder config containing the delegates to the strip measurement details.
sp – [in] Input space point used in the check.
spacepointPosition – [in] Spacepoint coordinates in xyz plane.
outputCoordinates – [out] The output vector to write to.
- Returns
Boolean that says if spacepoint is compatible with being inside the detector element.
-
template<typename ...R>
auto zip(R&&... r) Function that allows to zip some ranges to be used in a range-based for loop.
When wanting to mutate the entries, the result must be captured by value:
for(auto [a, b, c] : zip(ra, rb, rc)) { a+=2; }
Note
the behaviour is undefined if the ranges do not have equal range
- Template Parameters
R – The ranges type pack
- Parameters
r – The ranges parameter pack
Variables
-
static constexpr TrackIndexType kTrackIndexInvalid = std::numeric_limits<TrackIndexType>::max()
-
template<typename fitter>
constexpr bool LinearizerConcept = Acts::Concepts::Linearizer::LinearizerConcept<fitter>::value
-
constexpr int PolygonDynamic = -1
Tag to trigger specialization of a dynamic polygon.
-
static const std::vector<BinningValue> s_binningValues = {binX, binY, binZ, binR, binPhi, binRPhi, binH, binEta, binMag}
static list of all binning values
-
static constexpr ActsScalar s_curvilinearProjTolerance = 0.999995
Tolerance for not being within curvilinear projection this allows using the same curvilinear frame to eta = 6, validity tested with IntegrationTests/PropagationTest.
-
static constexpr ActsScalar s_epsilon = 3 * std::numeric_limits<ActsScalar>::epsilon()
Tolerance for being numerical equal for geometry building.
-
static const InfiniteBounds s_noBounds = {}
-
static constexpr ActsScalar s_onSurfaceTolerance = 1e-4
Tolerance for being on Surface.
Note
This is intentionally given w/o an explicit unit to avoid having to include the units header unnecessarily. With the native length unit of mm this corresponds to 0.1um.
-
static const Transform3 s_planeXY = Transform3::Identity()
-
static const Transform3 s_planeYZ = AngleAxis3(0.5 * M_PI, Vector3::UnitY()) * AngleAxis3(0.5 * M_PI, Vector3::UnitZ()) * Transform3::Identity()
-
static const Transform3 s_planeZX = AngleAxis3(-0.5 * M_PI, Vector3::UnitX()) * AngleAxis3(-0.5 * M_PI, Vector3::UnitZ()) * Transform3::Identity()
-
static ViewConfig s_viewFiltered = ViewConfig({255, 255, 0})
-
static ViewConfig s_viewGrid = ViewConfig({220, 0, 0})
-
static ViewConfig s_viewLine = ViewConfig({0, 0, 220})
-
static ViewConfig s_viewMeasurement = ViewConfig({255, 102, 0})
-
static ViewConfig s_viewParameter = ViewConfig({0, 0, 255})
-
static ViewConfig s_viewPassive = ViewConfig({240, 280, 0})
-
static ViewConfig s_viewPredicted = ViewConfig({51, 204, 51})
-
static ViewConfig s_viewSensitive = ViewConfig({0, 180, 240})
-
static ViewConfig s_viewSmoothed = ViewConfig({0, 102, 255})
-
static ViewConfig s_viewVolume = ViewConfig({220, 220, 0})
-
template<typename accessor>
constexpr bool SourceLinkAccessorConcept = Acts::Concepts::SourceLinkAccessor::SourceLinkAccessorConcept<accessor>::value
- template<typename stepper, typename state = typename stepper::State><
-
using ActsDynamicMatrix = Eigen::Matrix<ActsScalar, Eigen::Dynamic, Eigen::Dynamic>