ATLAS Offline Software
Loading...
Searching...
No Matches
MeasurementSelector.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4#pragma once
5
6// Alternative measurement selector
7//
8// This measurement selector is assuming the following
9// - the number of selected measurements is small (~<10)
10// - the number of measurement candidates is typically >> the
11// number of selected measuremebts
12// - the total number of candidate measurements can be large
13// - there is a simple one-to-one relation between bound state
14// parameters and measurement coordinates.
15
16#include "Acts/Utilities/Result.hpp"
17#include "Acts/Utilities/Delegate.hpp"
18#include "Acts/EventData/SourceLink.hpp"
19#include "Acts/TrackFinding/CombinatorialKalmanFilterError.hpp"
20#include "Acts/Surfaces/Surface.hpp"
21#include "Acts/Geometry/GeometryHierarchyMap.hpp"
22#include "Acts/EventData/Types.hpp"
23#include "Acts/EventData/TrackStatePropMask.hpp"
24#include "boost/container/small_vector.hpp"
25
26// for BaseTypes
27#include "Acts/EventData/TrackStateProxy.hpp"
28#include "Acts/Utilities/CalibrationContext.hpp"
29#include "Acts/EventData/TrackParameters.hpp"
30
31// for MeasurementSizeMax
32#include "Acts/EventData/MultiTrajectory.hpp"
33
34#include <utility>
35#include <type_traits>
36
37// Types to be used during measurement selection for the prediction and the
38// measurement for calibrated measurements after selection if the actual calibration is
39// executed after the calibration, and the trajectory, track state, bound parameter
40// types.
41template <typename derived_t>
43{
44 // the measurement type after the selection e.g. a Matrix<N,1>
45 template <std::size_t N>
46 using CalibratedMeasurement = typename Acts::detail_lt::FixedSizeTypes<N>::Coefficients;
47
48 // the measurement covariance type after the selection e.g. a Matrix<N,N>
49 template <std::size_t N>
50 using CalibratedMeasurementCovariance = typename Acts::detail_lt::FixedSizeTypes<N>::Covariance;
51
52 // the measurement type before the selection e.g. an Eigen::Map< Matrix<N,1> > if
53 // the calibration is performed after the selection
54 template <std::size_t N>
55 using PreSelectionMeasurement = typename Acts::detail_lt::FixedSizeTypes<N>::Coefficients;
56
57 // the measurement covariance type before the selection e.g. an Eigen::Map<Matrix<N,N> > if
58 // the calibration is performed after the selection
59 template <std::size_t N>
60 using PreSelectionMeasurementCovariance = typename Acts::detail_lt::FixedSizeTypes<N>::Covariance;
61
62 // e.g. the same as CalibratedMeasurement
63 template <std::size_t N>
64 using Predicted = typename Acts::detail_lt::FixedSizeTypes<N>::Coefficients;
65
66 // e.g. the same as CalibratedMeasurementCovariance
67 template <std::size_t N>
68 using PredictedCovariance = typename Acts::detail_lt::FixedSizeTypes<N>::Covariance;
69
70 // e.g. helper template to get the value_type from the container type
71 // e.g. helper template to get the value_type from the measurement range iterator type
72
73 template <typename T_MeasurementRangeIterator>
75 using value_type = typename T_MeasurementRangeIterator::value_type;
76 };
77
78 // abstract measurement range
79 // implements empty()
80 // and can be converted by the derived_t into an iterable range
81 // over measurements of concrete measurement container_type.
82 using abstract_measurement_range_t = std::ranges::iota_view<unsigned int, unsigned int>;
83
84 // the trajectory type to which states for selected measurements are to be added
85 using trajectory_t = typename derived_t::traj_t;
86 // the track state type for new track states
87 using TrackStateProxy = trajectory_t::TrackStateProxy;
88
89 // the value type usd for matrices
90 using MatrixFloatType = double;
91 using BoundTrackParameters = Acts::BoundTrackParameters;
92 using BoundMatrix = Acts::BoundMatrix;
93
94 using BoundState = std::tuple<BoundTrackParameters, BoundMatrix, double>;
95
96 // maximum dimension of measurements
97 static const std::size_t s_dimMax = 3;
98
99 // must be the same as what is used for the CKF
100 // @TODO how to get rid of this ?
101 static constexpr std::size_t s_maxBranchesPerSurface = 10;
102};
103
104//
108 template <class T_Matrix>
109 static constexpr std::size_t matrixColumns() { return T_Matrix::ColsAtCompileTime;}
110 template <class T_Matrix>
111 static constexpr std::size_t matrixRows() { return T_Matrix::RowsAtCompileTime;}
112
113 template <typename T_Float, class T_Matrix>
114 static auto matrixTypeCast(const T_Matrix &matrix) { return matrix.template cast<T_Float>(); }
115
116 template <class T_Matrix>
117 static auto transpose(const T_Matrix &matrix) { return matrix.transpose(); }
118
119 template <class T_Matrix>
120 static auto invert(const T_Matrix &matrix) { return matrix.inverse(); }
121};
122
123
124// Map from the measurement to bound state domain
125// it is assumed that there is a simple unambiguous association
126// between coordinates of the measurement domain and the bound state domain
127// e.g. measurement coordinate 1 maps to loc0 of the bound state.
128// @TODO use FixedSizeSubspace ? currently does not provide methods to directly
129// create the a sub-space covariance matrix from a full covariance matrix,
130// and to create the projector bitset, which is stored in the TrackState.
132
133 template <std::size_t N>
134 using type = std::array<unsigned char, N>;
135
136 template <std::size_t N>
137 static constexpr type<N> identity() {
138 type<N> ret;
139 for(int i=0; i<N; ++i) {
140 ret[i]=i;
141 }
142 return ret;
143 }
144};
145
146// utility to "project" a bound state parameter vector or covariance matrix onto the 1,2, ... N dimensional measurement domain
147// @TODO allow to influence resulting matrix type ?
148template <std::size_t N,class T_ResultType,class T_Matrix>
149T_ResultType project(ParameterMapping::type<N> parameter_map, const T_Matrix &matrix)
150{
151 using MatrixIndexMapType = unsigned char; // "char" to reduce the size of the map, and if not wide enough this entire
152 // concept is likely inefficient.
153 using MatrixIndexType = unsigned int; // @TODO or std::size_t ? does not matter
154
155 // ensure that index types are wide enough
156 static_assert( MeasurementSelectorMatrixTraits::matrixRows<T_Matrix>() < std::numeric_limits<MatrixIndexMapType>::max());
157 static_assert( N*MeasurementSelectorMatrixTraits::matrixRows<T_Matrix>() < std::numeric_limits<MatrixIndexType>::max());
158
159 T_ResultType ret;
161 // handle projection of paramteter vector
162 for (MatrixIndexType meas_i=0; meas_i<N; ++meas_i) {
163 assert( meas_i < parameter_map.size() );
164 ret(meas_i,0) = matrix( parameter_map[meas_i], 0);
165 }
166 }
167 else {
168 // handle projection of covariance matrix
169 // "project" matrix
170 for (MatrixIndexType meas_i=0; meas_i<N; ++meas_i) {
171 assert( meas_i < parameter_map.size());
172 MatrixIndexType param_i = parameter_map[meas_i];
173 for (MatrixIndexType meas_j=0; meas_j<N; ++meas_j) {
174 assert( meas_j < parameter_map.size());
175 ret(meas_i,meas_j) = matrix(param_i, parameter_map[meas_j]);
176 }
177 }
178 }
179 return ret;
180}
181
182
183
184
185// helper method to compute a chi2 for the difference of two "measurement" and covariance pairs
186template <typename measurement_vector_t, typename measurement_cov_matrix_t,
187 typename predicted_vector_t, typename predicted_cov_matrix_t>
188double computeChi2(const measurement_vector_t &a,
189 const measurement_cov_matrix_t &a_cov,
190 const predicted_vector_t &b,
191 const predicted_cov_matrix_t &b_cov) {
192
193 // just sum sanity checks that a and b have the correct dimensions and that
194 // the chi2 can actually be computed.
196 static_assert( MeasurementSelectorMatrixTraits::matrixColumns<predicted_vector_t>() == 1); // b is vector
207
208 // @TODO remove abstraction i.e. assume matrix has the interface of an eigen matrix
209 auto inv_ab_cov( MeasurementSelectorMatrixTraits::invert(a_cov+b_cov) );
210 auto diff( a-b);
211 return (MeasurementSelectorMatrixTraits::transpose(diff) * inv_ab_cov * diff)(0,0);
212}
213
214// Collection to hold the n-"best" candidates
215// The objects of type PayloadType must support assignment operation, and must
216// be default constructible. Moreover it must be possible to provide a
217// "comparison" operator to order the payload objects.
218template <std::size_t N, class PayloadType >
220 using IndexType = unsigned short; // @TODO or char ? If N>>10 this concept is likely
221 // inefficient
222 // using PayloadType = Payload<DIM>;
223 TopCollection(std::size_t max_n) {
224 init(max_n);
225 }
226
227 // @param max_n the maximum number of top-candidates is fixed by the template parameter
228 // N but can be reduced further to this number
229 void init(std::size_t max_n) {
230 assert( max_n < N);
231 m_nextSlot=0;
232 m_maxSlots=max_n;
233 m_order[0]=0;
234 }
235 // @param get a slot to hold a new candidate which is not necessarily accepted in the list
236 // of the n-top candidates
237 PayloadType &slot() {
238 return m_slots[m_order[m_nextSlot] ];
239 }
240 // @param idx get the specified filled slot (read only) indicated by the index, where the index does not
241 // indicate the order in the top-candidate list
242 const PayloadType &getSlot(IndexType idx) const {
243 return m_slots[idx];
244 }
245 // @param idx get the specified filled slot indicated by the index, where the index does not
246 // indicate the order in the top-candidate list
247 PayloadType &getSlot(IndexType idx) {
248 return m_slots[idx];
249 }
250 // @param test whether the given index points to one of the accepted top candidates
251 bool isValid(IndexType idx) const {
252 return idx < m_nextSlot;
253 }
254 // Accept the element of the latest slot provided there is still a free slot or it is better than
255 // the worst element.
256 // @param comparison operator to compute e.g. a<b where "smaller" means "better"
257 void acceptAndSort(std::function<bool(const PayloadType &a, const PayloadType &b)> comparison) {
258 // bubble best element to top
259 for (unsigned int slot_i = m_nextSlot;
260 slot_i-- > 0
261 && !comparison(m_slots[ m_order[slot_i] ],m_slots[ m_order[slot_i+1] ]);) {
262 std::swap(m_order[slot_i],m_order[slot_i+1]);
263 }
264 // if there are still free slot increase the number of used slots
265 if (m_nextSlot < m_maxSlots) {
266 ++m_nextSlot;
268 }
269 }
270 // Accept the element of the latest slot provided there is still a free slot.
272 // if there are still free slot increase the number of used slots
273 if (m_nextSlot < m_maxSlots) {
274 ++m_nextSlot;
276 }
277 return (m_nextSlot < m_maxSlots);
278 }
279
280 bool empty() const {
281 return m_nextSlot==0;
282 }
283
284 IndexType size() const { return m_nextSlot; }
285
286 // helper to iterate over the slot indices in sorting order from best to worst
287 typename std::array<IndexType, N+1>::const_iterator begin() { return m_order.begin(); }
288 typename std::array<IndexType, N+1>::const_iterator end() { return m_order.begin()+m_nextSlot; }
289
290 std::array<PayloadType, N+1> m_slots; // storage for the slots
291 std::array<IndexType, N+1> m_order; // order of the filled slots
292 IndexType m_nextSlot = 0; // the index of the next free slot
293 IndexType m_maxSlots = 0; // maximum number of top-slots
294};
295
296
297
300 std::vector<float> etaBins{};
302 std::vector<std::pair<float, float> > chi2CutOff{ {15,25} };
304 std::vector<std::size_t> numMeasurementsCutOff{1};
305};
306
307// Measurement type specific measirement selector
308// Assumptions:
309// - begin and end source_link_iterator point to a contiguous range of measurements in a single container
310// - all measurements in this range have the same dimensionality.
311// - the mapping of the bound parameters to the measurement domain for each measurement in this range
312// is identical
313template <std::size_t NMeasMax,
314 std::size_t DIMMAX,
315 typename derived_t >
317
318 using Config = Acts::GeometryHierarchyMap<AtlasMeasurementSelectorCuts>;
320
322
323protected:
328
329 const derived_t &derived() const { return *static_cast<const derived_t *>(this); }
330
331 // helper to create a projector bitset from a map from bound parameters to coordinates
333 template <std::size_t N>
334 static
335 Acts::ProjectorBitset create(const ParameterMapping::type<N> &parameter_map) {
336 constexpr std::size_t nrows = Acts::MultiTrajectoryTraits::MeasurementSizeMax;
337 constexpr std::size_t ncols = Acts::eBoundSize;
338
339 std::bitset<nrows * ncols> proj_bitset {};
340
341 for (unsigned int col_i=0; col_i<N; ++col_i) {
342 unsigned int row_i = parameter_map[col_i];
343 unsigned int idx = col_i *nrows + row_i; // @TODO handle row major and column major correctly
344 proj_bitset[ (nrows * ncols - 1) - idx ] = 1;
345 }
346 return proj_bitset.to_ullong();
347 }
348 };
349
350 // helper to create states without the calibrated measurement, covariance, and sourcelink
351 // if more than one state is created the states after the first will share
352 // the jacobi and the prediction, if outlier_states is true no storage is created
353 // for filtered.
354 // @TODO should there be the possibility to have multiple states which mix outlier and non-outlier states ?
355 static void createStates(std::size_t n_new_states,
356 const T_BoundState& boundState,
357 std::size_t prevTip,
358 trajectory_t& trajectory,
359 const Acts::BoundSubspaceIndices& subspaceIndices,
360 boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface> &track_states,
361 const Acts::Logger& logger,
362 bool outlier_states) {
363 // create track states for all compatible measurement candidates
364 const auto &boundParams = derived_t::boundParams(boundState);
365 const auto &pathLength = derived_t::pathLength(boundState);
366
367 using PM = Acts::TrackStatePropMask;
368
369 track_states.reserve( n_new_states );
370 std::optional<TrackStateProxy> firstTrackState{
371 std::nullopt};
372 for (unsigned int state_i=0; state_i<n_new_states; ++state_i) {
373 PM mask = PM::Predicted | PM::Filtered | PM::Jacobian | PM::Calibrated;
374
375 if (firstTrackState.has_value()) {
376 // subsequent track states don't need storage for these as they will
377 // be shared
378 mask &= ~PM::Predicted & ~PM::Jacobian;
379 }
380
381 if (outlier_states) {
382 // outlier won't have separate filtered parameters
383 mask &= ~PM::Filtered;
384 }
385
386 TrackStateProxy trackState = trajectory.makeTrackState(mask, prevTip);
387 ACTS_VERBOSE("Create SourceLink output track state #"
388 << trackState.index() << " with mask: " << mask);
389
390 if (firstTrackState.has_value()) {
391 trackState.shareFrom(*firstTrackState, PM::Predicted);
392 trackState.shareFrom(*firstTrackState, PM::Jacobian);
393 }
394 else {
395 // only set these for first
396 trackState.predicted() = boundParams.parameters();
397 if (boundParams.covariance()) {
398 trackState.predictedCovariance() = *boundParams.covariance();
399 }
400 trackState.jacobian() = derived_t::boundJacobiMatrix(boundState);
401 firstTrackState = trackState;
402 }
403 trackState.pathLength() = pathLength;
404
405 trackState.setReferenceSurface(boundParams.referenceSurface().getSharedPtr());
406
407 trackState.setProjectorSubspaceIndices(subspaceIndices);
408
409 Acts::TrackStateType typeFlags = trackState.typeFlags();
410 if (trackState.referenceSurface().surfaceMaterial() != nullptr) {
411 typeFlags.set(Acts::TrackStateFlag::MaterialFlag);
412 }
413 typeFlags.set(Acts::TrackStateFlag::ParameterFlag);
414
415 // @TODO these track states still need some additional processing. Should there be a special
416 // flag for this ?
417 track_states.push_back( trackState.index());
418 }
419 }
420
421
422 // utility struct to temporarily store data about the best measurements
423 template <std::size_t DIM, typename T_SourceLink>
425 using MeasCovPair = std::pair < typename MeasurementSelectorTraits<derived_t>::template PreSelectionMeasurement<DIM>,
426 typename MeasurementSelectorTraits<derived_t>::template PreSelectionMeasurementCovariance<DIM> >;
428 std::optional<T_SourceLink> m_sourceLink;
429 float m_chi2;
431 };
432
433 // type and dimension specific function to select measurements from the range defined by the source link iterators.
434 // will iterate over the contiguous measurement range defined by the source link iterators where the measurements
435 // are contained in the given container. The selection lopp will get the measurement and covariance with the
436 // help of a preCalibrator. select measurements based on smallest chi2 wrt. the prediction, then optionally
437 // apply a full calibrator after the selection and finally create track states.
438 template <std::size_t DIM, typename T_MeasurementRange>
439 Acts::Result<boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface> >
440 selectMeasurementsCreateTrackStates(const Acts::GeometryContext& geometryContext,
441 const Acts::CalibrationContext& calibrationContext,
442 const Acts::Surface& surface,
443 const T_BoundState& boundState,
444 T_MeasurementRange &&measurement_range,
445 std::size_t prevTip,
446 trajectory_t& trajectory,
447 const Acts::Logger& logger,
448 const std::size_t numMeasurementsCut,
449 const std::pair<float,float>& maxChi2Cut,
450 bool forced) const {
451 Acts::Result<boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface> >
452 result = boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface>{};
453
454 using iterator_t = decltype(measurement_range.begin());
455 using container_value_t = typename MeasurementSelectorTraits<derived_t>::template MeasurementContainerTraits<iterator_t>::value_type;
456 using BaseElementType = std::remove_cv_t<std::remove_pointer_t< container_value_t > >;
457 // get calibrator
458 using TheMatchingMeasurement = MatchingMeasurement<DIM, container_value_t >;
459 using MeasCovPair = TheMatchingMeasurement::MeasCovPair;
460
461 using Predicted = typename MeasurementSelectorTraits<derived_t>::template Predicted<DIM>;
462 using PredictedCovariance = typename MeasurementSelectorTraits<derived_t>::template PredictedCovariance<DIM>;
463
464
465
466 // get prediction in the measurement domain
467 Acts::SubspaceIndices<DIM> parameter_map = derived().template parameterMap<DIM>(geometryContext,
468 calibrationContext,
469 surface);
470 auto predicted
471 = std::make_pair( project<DIM, Predicted>(parameter_map,
472 derived().boundParams(boundState).parameters()),
474 derived().boundParams(boundState).covariance().value()));
475
476 // select n measurents with the smallest chi2.
477 auto preCalibrator = derived().template preCalibrator<DIM, BaseElementType>();
478 auto postCalibrator = derived().template postCalibrator<DIM, BaseElementType>();
479 TopCollection<NMeasMax, TheMatchingMeasurement > selected_measurements(numMeasurementsCut);
480 for ( const auto &measurement : measurement_range ) {
481 TheMatchingMeasurement &matching_measurement=selected_measurements.slot();
482 if (forced && postCalibrator) {
483 // skip preCalibrator and computeChi2, by doing everything else. We start with what's done if no preCalibrator.
484 const auto &m = derived().forwardToCalibrator(measurement);
485 matching_measurement.m_measurement = std::make_pair( m.template localPosition<DIM>(), m.template localCovariance<DIM>() );
486 matching_measurement.m_chi2 = 0.0f;
487 matching_measurement.m_sourceLink=measurement;
488 if (!selected_measurements.acceptNoSort()) break;
489 } else {
490 matching_measurement.m_measurement = preCalibrator(geometryContext,
491 calibrationContext,
492 derived().forwardToCalibrator(measurement),
493 derived().boundParams(boundState));
494 matching_measurement.m_chi2 = computeChi2(matching_measurement.m_measurement.first,
495 matching_measurement.m_measurement.second,
496 predicted.first,
497 predicted.second);
498 // only consider measurements which pass the outlier chi2 cut
499 if (matching_measurement.m_chi2<maxChi2Cut.second) {
500 matching_measurement.m_sourceLink=measurement;
501 selected_measurements.acceptAndSort([](const TheMatchingMeasurement &a,
502 const TheMatchingMeasurement &b) {
503 return a.m_chi2 < b.m_chi2;
504 });
505 }
506 }
507 }
508
509 // apply final calibration to n-best measurements
510 using post_calib_meas_cov_pair_t
511 = std::pair<typename MeasurementSelectorTraits<derived_t>::template CalibratedMeasurement<DIM>,
512 typename MeasurementSelectorTraits<derived_t>::template CalibratedMeasurementCovariance<DIM> >;
513
514 // need extra temporary buffer if the types to store "measurements" during measurement selection
515 // and calibrated measurements after the measurement selection are not identical
516 static constexpr bool pre_and_post_calib_types_agree
517 = std::is_same< typename MeasurementSelectorTraits<derived_t>::template CalibratedMeasurement<DIM>,
518 typename MeasurementSelectorTraits<derived_t>::template PreSelectionMeasurement<DIM> >::value
519 && std::is_same< typename MeasurementSelectorTraits<derived_t>::template CalibratedMeasurementCovariance<DIM>,
520 typename MeasurementSelectorTraits<derived_t>::template PreSelectionMeasurementCovariance<DIM> >::value;
521
522 using Empty = struct {};
523 typename std::conditional< !pre_and_post_calib_types_agree,
524 std::array< post_calib_meas_cov_pair_t, NMeasMax>, // @TODO could use boost static_vector instead
525 Empty>::type calibrated;
526 if (postCalibrator) {
527 typename std::conditional<!pre_and_post_calib_types_agree, unsigned int, Empty>::type calibrated_meas_cov_i{};
528
529 // apply final calibration, recompute chi2
531 idx: selected_measurements) {
532 TheMatchingMeasurement &a_selected_measurement = selected_measurements.getSlot(idx);
533
534 // helper to select the destination storage which is the extra temporary buffer
535 // if the measurement storage types during and after selection are different.
536 post_calib_meas_cov_pair_t &calibrated_measurement
537 = [&calibrated, &a_selected_measurement, calibrated_meas_cov_i]() -> post_calib_meas_cov_pair_t & {
538 if constexpr(pre_and_post_calib_types_agree) {
539 (void) calibrated;
540 (void) calibrated_meas_cov_i;
541 return a_selected_measurement.m_measurement;
542 }
543 else {
544 (void) a_selected_measurement;
545 assert(calibrated_meas_cov_i < calibrated.size());
546 return calibrated[calibrated_meas_cov_i];
547 }
548 }();
549
550 // apply the calibration
551 calibrated_measurement = postCalibrator(geometryContext,
552 calibrationContext,
553 derived().forwardToCalibrator(a_selected_measurement.m_sourceLink.value()),
554 derived().boundParams(boundState));
555 // update chi2 using calibrated measurement
556 a_selected_measurement.m_chi2 = computeChi2(calibrated_measurement.first,
557 calibrated_measurement.second,
558 predicted.first,
559 predicted.second);
560 // ... and set outlier flag
561 a_selected_measurement.m_isOutLier = (!forced && a_selected_measurement.m_chi2 >= maxChi2Cut.first);
562 if constexpr(!pre_and_post_calib_types_agree) {
563 ++calibrated_meas_cov_i;
564 }
565 }
566 }
567 else {
568 // if no final calibration is performed only the outlier flag still needs to be set
570 idx: selected_measurements) {
571 TheMatchingMeasurement &a_selected_measurement = selected_measurements.getSlot(idx);
572 a_selected_measurement.m_isOutLier = (!forced && a_selected_measurement.m_chi2 >= maxChi2Cut.first);
573 }
574 }
575
576 // First Create states without setting information about the calibrated measurement for the selected measurements
577 // @TODO first create state then copy measurements, or crete state by state and set measurements ?
578 // the lastter has the "advantage" that the outlier flag can be set individually
579 // the former has the advantage that part of the state creation code is independent of the
580 // the measurement.
581
582 Acts::BoundSubspaceIndices boundSubspaceIndices;
583 std::copy(parameter_map.begin(), parameter_map.end(), boundSubspaceIndices.begin());
584 createStates( selected_measurements.size(),
585 boundState,
586 prevTip,
587 trajectory,
588 boundSubspaceIndices,
589 *result,
590 logger,
591 (!selected_measurements.empty()
592 ? selected_measurements.getSlot( *(selected_measurements.begin())).m_isOutLier
593 : false) );
594 assert( result->size() == selected_measurements.size() );
595
596 // helper to determine whether calibrated storeage is to be used
597 auto use_calibrated_storage = [&postCalibrator]() -> bool {
598 if constexpr(pre_and_post_calib_types_agree) {
599 (void) postCalibrator;
600 return false;
601 }
602 else {
603 // this is only known during runtime
604 // but it should not be tested if the types agree.
605 return postCalibrator;
606 }
607 };
608
609 // copy selected measurements to pre-created states
610 unsigned int state_i=0;
612 idx: selected_measurements) {
613 assert( state_i < result->size());
614 TrackStateProxy trackState( trajectory.getTrackState( (*result)[state_i] ) );
615 TheMatchingMeasurement &a_selected_measurement = selected_measurements.getSlot(idx);
616 trackState.setUncalibratedSourceLink(derived().makeSourceLink(std::move(a_selected_measurement.m_sourceLink.value())));
617 // flag outliers accordingly, so that they are handled correctly by the post processing
618 trackState.typeFlags().set( a_selected_measurement.m_isOutLier
619 ? Acts::TrackStateFlag::OutlierFlag
620 : Acts::TrackStateFlag::MeasurementFlag );
621 trackState.allocateCalibrated(DIM);
622 if (use_calibrated_storage()) {
623 // if the final clibration is performed after the selection then
624 // copy these measurements and covariances to the track states
625 assert( use_calibrated_storage() == !pre_and_post_calib_types_agree);
626 if constexpr(!pre_and_post_calib_types_agree) {
627 assert( state_i < calibrated.size());
628 trackState.template calibrated<DIM>()
630 trackState.template calibratedCovariance<DIM>()
632 trackState.chi2() = a_selected_measurement.m_chi2;
633 }
634 }
635 else {
636 trackState.template calibrated<DIM>()
638 trackState.template calibratedCovariance<DIM>()
640 trackState.chi2() = a_selected_measurement.m_chi2;
641 }
642 ++state_i;
643 }
644 return result;
645 }
646
647 template <typename parameters_t>
648 static std::size_t getEtaBin(const parameters_t& boundParameters,
649 const std::vector<float> &etaBins) {
650 if (etaBins.empty()) {
651 return 0u; // shortcut if no etaBins
652 }
653 const float eta = std::abs(std::atanh(std::cos(boundParameters.parameters()[Acts::eBoundTheta])));
654 std::size_t bin = 0;
655 for (auto etaBin : etaBins) {
656 if (etaBin >= eta) {
657 break;
658 }
659 bin++;
660 }
661 return bin;
662 }
663
664public:
665 // get numMeasurement and maxChi2 cuts for the given surface
666 // @TODO should the cuts just depend on the surface or really on the bound parameters (i.e. bound eta)
667 std::tuple< std::size_t, std::pair<float,float> > getCuts(const Acts::Surface& surface,
668 const T_BoundState& boundState,
669 const Acts::Logger& logger) const {
670 std::tuple< std::size_t, std::pair<float,float> > result;
671 std::size_t &numMeasurementsCut = std::get<0>(result);
672 std::pair<float,float> &maxChi2Cut = std::get<1>(result);
673 // Get geoID of this surface
674 auto geoID = surface.geometryId();
675 // Find the appropriate cuts
676 auto cuts = m_config.find(geoID);
677 if (cuts == m_config.end()) {
678 // indicats failure
679 numMeasurementsCut = 0;
680 }
681 else {
682 // num measurement Cut
683 // getchi2 cut
684 std::size_t eta_bin = getEtaBin(derived().boundParams(boundState), cuts->etaBins);
685 numMeasurementsCut = (!cuts->numMeasurementsCutOff.empty()
686 ? cuts->numMeasurementsCutOff[ std::min(cuts->numMeasurementsCutOff.size()-1, eta_bin) ]
687 : NMeasMax);
688 maxChi2Cut = ! cuts->chi2CutOff.empty()
689 ? cuts->chi2CutOff[ std::min(cuts->chi2CutOff.size()-1, eta_bin) ]
690 : std::make_pair<float,float>(std::numeric_limits<float>::max(),
691 std::numeric_limits<float>::max());
692 ACTS_VERBOSE("Get cut for eta-bin="
693 << (eta_bin < cuts->etaBins.size() ? cuts->etaBins[eta_bin] : std::numeric_limits<float>::max())
694 << ": chi2 (max,max-outlier) " << maxChi2Cut.first << ", " << maxChi2Cut.second
695 << " max.meas." << numMeasurementsCut);
696 }
697 return result;
698 }
699};
700
701// Measurement selector which calls a type specific selection method
702// the measurement selector expects a list of measurement containers to be provided
703// by the source link iterators as well as the index of the container to be used for the
704// given measurement range, which is defined by the source link iterators.
705// the measurement range must refer to a contiguous range of measurements contained
706// within a single measurement container.
707// The "container" itself is a variant of particular measurement containers with associated
708// dimension i.e. number of coordinates per measurement. A member funcion specific to one of
709// the alternatives of the variant will be called to perform the selection and track state
710// creation.
711template <std::size_t NMeasMax,
712 typename derived_t,
713 typename measurement_container_variant_t >
714struct MeasurementSelectorWithDispatch : public MeasurementSelectorBase< NMeasMax, MeasurementSelectorTraits<derived_t>::s_dimMax, derived_t> {
715
717 using T_BoundState = typename Base::T_BoundState;
718 using trajectory_t = typename Base::trajectory_t;
719 using TrackStateProxy = typename Base::TrackStateProxy;
720 static constexpr std::size_t s_maxBranchesPerSurface = Base::s_maxBranchesPerSurface;
721
722 // helper to get the maximum number of measurement diemsions.
723 static constexpr std::size_t dimMax() {
725 }
726
727 Acts::Result<boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface> >
728 createTrackStates(const Acts::GeometryContext& geometryContext,
729 const Acts::CalibrationContext& calibrationContext,
730 const Acts::Surface& surface,
731 const T_BoundState& boundState,
732 typename TrackStateProxy::IndexType prevTip,
733 [[maybe_unused]] std::vector<TrackStateProxy>& trackStateCandidates,
734 trajectory_t& trajectory,
735 const Acts::Logger& logger) const {
736
737 Acts::Result<boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface> >
738 result = Acts::CombinatorialKalmanFilterError::MeasurementSelectionFailed;
739 // get associated measurement container and the relevant measurement range for the given surface.
740 auto [a_measurement_container_variant_ptr, range, forced] = this->derived().containerAndRange(surface);
741 if (!forced && !this->derived().expectMeasurements( surface, a_measurement_container_variant_ptr, range)) {
742 result = result.failure(Acts::CombinatorialKalmanFilterError::NoMeasurementExpected);
743 return result;
744 }
745 if (!range.empty()) {
746 auto [numMeasurementsCut, maxChi2Cut] = this->getCuts(surface,boundState, logger);
747 // numMeasurementsCut is == 0 in case getCuts failed
748 // anyway cannot select anything if numMeasurementsCut==0;
749 if (numMeasurementsCut>0) {
750
751 result = std::visit( [this,
752 &geometryContext,
753 &calibrationContext,
754 &surface,
755 &boundState,
756 &range,
757 prevTip,
758 &trajectory,
759 &logger,
760 numMeasurementsCut,
761 &maxChi2Cut,
762 forced] (const auto &measurement_container_with_dimension) {
763 using ArgType = std::remove_cv_t<std::remove_reference_t< decltype(measurement_container_with_dimension) > >;
764 constexpr std::size_t DIM = ArgType::dimension();
765 auto measurement_range = this->derived().rangeForContainer(measurement_container_with_dimension,range);
766 return
767 this->template
769 calibrationContext,
770 surface,
771 boundState,
772 std::move(measurement_range),
773 prevTip,
774 trajectory,
775 logger,
776 numMeasurementsCut,
777 maxChi2Cut,
778 forced);
779 },
780 *a_measurement_container_variant_ptr);
781 }
782 }
783 else {
784 result = Acts::Result<boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface> >::success({});
785 }
786 return result;
787 }
788};
789
790template <std::size_t NMeasMax,
791 typename derived_t,
792 typename measurement_container_variant_t>
793struct MeasurementSelectorBaseImpl : public MeasurementSelectorWithDispatch<NMeasMax, derived_t, measurement_container_variant_t> {
796
797 // get the bound parameters and covariance of the bound state
799 return std::get<0>(boundState);
800 }
801
802 // get the jacobi matrix of the bound state
804 return std::get<1>(boundState);
805 }
806 // get the accumulated path length of the bound state
807 static double pathLength(const T_BoundState &boundState) {
808 return std::get<2>(boundState);
809 }
810
811 // create a source link from the measurement
812 template <typename T_Value>
813 static Acts::SourceLink makeSourceLink(T_Value &&value) {
814 return Acts::SourceLink{value};
815 }
816
817 // perform simple transformation to cerate the type
818 // that is passed to the calibrator.
819 // by default pass the measurement by reference not pointer
820 template <typename T>
821 static const auto &forwardToCalibrator(const T &a) {
822 if constexpr( std::is_same<T, std::remove_pointer_t<T> >::value ) {
823 return a;
824 }
825 else {
826 return *a;
827 }
828 }
829
830
831 // get mapping between bound state parameters and coordinates in the measurement domain
832 // in most cases this is just the identity operation e.g. loc0 -> coord0 and loc1 -> coord1
833 // i.e. map[0]=0; map[1]=1;
834 // @TODO should this be measurement type specific, or is dimension and the surface good enough ?
835 template <std::size_t DIM>
837 parameterMap(const Acts::GeometryContext&,
838 const Acts::CalibrationContext&,
839 const Acts::Surface&) {
841 }
842
843 // By default the methods postCalibrator and preCalibrator return delegates:
844 template <std::size_t DIM, typename measurement_t>
845 using PreCalibrator = Acts::Delegate<
846 std::pair < typename MeasurementSelectorTraits<derived_t>::template PreSelectionMeasurement<DIM>,
847 typename MeasurementSelectorTraits<derived_t>::template PreSelectionMeasurementCovariance<DIM> >
848 (const Acts::GeometryContext&,
849 const Acts::CalibrationContext&,
850 const measurement_t &,
852
853 // Since the measurement types used during measurement selection and after measurement selection
854 // might be different so are the types of the calibrator delegates
855 template <std::size_t DIM, typename measurement_t>
856 using PostCalibrator = Acts::Delegate<
857 std::pair < typename MeasurementSelectorTraits<derived_t>::template CalibratedMeasurement<DIM>,
858 typename MeasurementSelectorTraits<derived_t>::template CalibratedMeasurementCovariance<DIM> >
859 (const Acts::GeometryContext&,
860 const Acts::CalibrationContext&,
861 const measurement_t &,
863
867 template <std::size_t DIM, typename measurement_t>
869 postCalibrator() const; // not implemented
870
871 // the "calibrator" which is used during the measuremnt selection
873 // this delegate must be connected or be a valid functor. If not this likely will lead to an
874 // exception or seg fault.
875 template <std::size_t DIM, typename measurement_t>
877 preCalibrator() const; // not implemented
878
879
884 std::tuple<const measurement_container_variant_t &, abstract_measurement_range_t>
885 containerAndRange(const Acts::Surface &surface) const; // not implemented
886
891 bool expectMeasurements([[maybe_unused]] const Acts::Surface &surface,
892 [[maybe_unused]] const measurement_container_variant_t *container_variant_ptr,
893 const abstract_measurement_range_t &abstract_range) const;
894
899 template <typename measurement_container_t>
900 static auto
901 rangeForContainer(const measurement_container_t &concrete_container, const abstract_measurement_range_t &abstract_range);
902};
Scalar eta() const
pseudorapidity method
static Double_t a
double computeChi2(const measurement_vector_t &a, const measurement_cov_matrix_t &a_cov, const predicted_vector_t &b, const predicted_cov_matrix_t &b_cov)
T_ResultType project(ParameterMapping::type< N > parameter_map, const T_Matrix &matrix)
void diff(const Jet &rJet1, const Jet &rJet2, std::map< std::string, double > varDiff)
Difference between jets - Non-Class function required by trigger.
Definition Jet.cxx:631
static Root::TMsgLogger logger("iLumiCalc")
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)
std::vector< std::size_t > numMeasurementsCutOff
Maximum number of associated measurements on a single surface.
std::vector< std::pair< float, float > > chi2CutOff
Maximum local chi2 contribution.
std::vector< float > etaBins
bins in |eta| to specify variable selections
MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::abstract_measurement_range_t abstract_measurement_range_t
static const MeasurementSelectorTraits< derived_t >::BoundMatrix & boundJacobiMatrix(const T_BoundState &boundState)
MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::BoundState T_BoundState
bool expectMeasurements(const Acts::Surface &surface, const measurement_container_variant_t *container_variant_ptr, const abstract_measurement_range_t &abstract_range) const
const PreCalibrator< DIM, measurement_t > & preCalibrator() const
const PostCalibrator< DIM, measurement_t > & postCalibrator() const
the calibrator used after the measurement selection which does not have to be "connected"
static double pathLength(const T_BoundState &boundState)
std::tuple< const measurement_container_variant_t &, abstract_measurement_range_t > containerAndRange(const Acts::Surface &surface) const
Get container and index range for measurements on the given surface.
Acts::Delegate< std::pair< typename MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::template PreSelectionMeasurement< DIM >, typename MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::template PreSelectionMeasurementCovariance< DIM > >(const Acts::GeometryContext &, const Acts::CalibrationContext &, const measurement_t &, const typename MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::BoundTrackParameters &)> PreCalibrator
static const auto & forwardToCalibrator(const T &a)
Acts::Delegate< std::pair< typename MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::template CalibratedMeasurement< DIM >, typename MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::template CalibratedMeasurementCovariance< DIM > >(const Acts::GeometryContext &, const Acts::CalibrationContext &, const measurement_t &, const typename MeasurementSelectorTraits< AtlasMeasurementSelector< NMeasMax, traj_t, measurement_container_variant_t > >::BoundTrackParameters &)> PostCalibrator
static Acts::SourceLink makeSourceLink(T_Value &&value)
static const MeasurementSelectorTraits< derived_t >::BoundTrackParameters & boundParams(const T_BoundState &boundState)
static auto rangeForContainer(const measurement_container_t &concrete_container, const abstract_measurement_range_t &abstract_range)
Create a range over elements of the given container from an abstract range.
ParameterMapping::type< DIM > parameterMap(const Acts::GeometryContext &, const Acts::CalibrationContext &, const Acts::Surface &)
std::pair< typename MeasurementSelectorTraits< derived_t >::template PreSelectionMeasurement< DIM >, typename MeasurementSelectorTraits< derived_t >::template PreSelectionMeasurementCovariance< DIM > > MeasCovPair
static Acts::ProjectorBitset create(const ParameterMapping::type< N > &parameter_map)
const derived_t & derived() const
MeasurementSelectorTraits< derived_t > traits
std::tuple< std::size_t, std::pair< float, float > > getCuts(const Acts::Surface &surface, const T_BoundState &boundState, const Acts::Logger &logger) const
static void createStates(std::size_t n_new_states, const T_BoundState &boundState, std::size_t prevTip, trajectory_t &trajectory, const Acts::BoundSubspaceIndices &subspaceIndices, boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface > &track_states, const Acts::Logger &logger, bool outlier_states)
static std::size_t getEtaBin(const parameters_t &boundParameters, const std::vector< float > &etaBins)
typename MeasurementSelectorTraits< derived_t >::trajectory_t trajectory_t
Acts::Result< boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface > > selectMeasurementsCreateTrackStates(const Acts::GeometryContext &geometryContext, const Acts::CalibrationContext &calibrationContext, const Acts::Surface &surface, const T_BoundState &boundState, T_MeasurementRange &&measurement_range, std::size_t prevTip, trajectory_t &trajectory, const Acts::Logger &logger, const std::size_t numMeasurementsCut, const std::pair< float, float > &maxChi2Cut, bool forced) const
typename MeasurementSelectorTraits< derived_t >::BoundState T_BoundState
typename MeasurementSelectorTraits< derived_t >::TrackStateProxy TrackStateProxy
Acts::GeometryHierarchyMap< AtlasMeasurementSelectorCuts > Config
static constexpr std::size_t matrixRows()
static constexpr std::size_t matrixColumns()
matrix adapter for Eigen additionally need +,- and *
static auto transpose(const T_Matrix &matrix)
static auto matrixTypeCast(const T_Matrix &matrix)
static auto invert(const T_Matrix &matrix)
typename T_MeasurementRangeIterator::value_type value_type
Acts::BoundTrackParameters BoundTrackParameters
trajectory_t::TrackStateProxy TrackStateProxy
typename Acts::detail_lt::FixedSizeTypes< N >::Coefficients PreSelectionMeasurement
std::ranges::iota_view< unsigned int, unsigned int > abstract_measurement_range_t
typename Acts::detail_lt::FixedSizeTypes< N >::Covariance PreSelectionMeasurementCovariance
typename Acts::detail_lt::FixedSizeTypes< N >::Covariance PredictedCovariance
typename Acts::detail_lt::FixedSizeTypes< N >::Coefficients Predicted
static constexpr std::size_t s_maxBranchesPerSurface
typename derived_t::traj_t trajectory_t
static const std::size_t s_dimMax
typename Acts::detail_lt::FixedSizeTypes< N >::Covariance CalibratedMeasurementCovariance
std::tuple< BoundTrackParameters, BoundMatrix, double > BoundState
typename Acts::detail_lt::FixedSizeTypes< N >::Coefficients CalibratedMeasurement
Acts::Result< boost::container::small_vector< typename TrackStateProxy::IndexType, s_maxBranchesPerSurface > > createTrackStates(const Acts::GeometryContext &geometryContext, const Acts::CalibrationContext &calibrationContext, const Acts::Surface &surface, const T_BoundState &boundState, typename TrackStateProxy::IndexType prevTip, std::vector< TrackStateProxy > &trackStateCandidates, trajectory_t &trajectory, const Acts::Logger &logger) const
typename Base::trajectory_t trajectory_t
MeasurementSelectorBase< NMeasMax, MeasurementSelectorTraits< derived_t >::s_dimMax, derived_t > Base
static constexpr std::size_t dimMax()
typename Base::T_BoundState T_BoundState
static constexpr std::size_t s_maxBranchesPerSurface
typename Base::TrackStateProxy TrackStateProxy
static constexpr type< N > identity()
std::array< unsigned char, N > type
std::array< PayloadType, N+1 > m_slots
std::array< IndexType, N+1 >::const_iterator begin()
std::array< IndexType, N+1 > m_order
const PayloadType & getSlot(IndexType idx) const
unsigned short IndexType
IndexType size() const
PayloadType & slot()
TopCollection(std::size_t max_n)
bool isValid(IndexType idx) const
void acceptAndSort(std::function< bool(const PayloadType &a, const PayloadType &b)> comparison)
std::array< IndexType, N+1 >::const_iterator end()
PayloadType & getSlot(IndexType idx)
void init(std::size_t max_n)