ATLAS Offline Software
Loading...
Searching...
No Matches
TwoStateCombiner.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
11
16//
18
19namespace {
20 /*
21 * Utilities used for the forward / backward
22 * smoother combination.
23 */
24void
25combineWithWeight(Trk::ComponentParameters& mergeTo,
26const Trk::ComponentParameters& addThis)
27{
28 const Trk::TrackParameters* firstTrackParameters = mergeTo.params.get();
29 const AmgVector(5)& firstParameters = firstTrackParameters->parameters();
30 double const firstWeight = mergeTo.weight;
31
32 const Trk::TrackParameters* secondTrackParameters = addThis.params.get();
33 const AmgVector(5)& secondParameters = secondTrackParameters->parameters();
34 double const secondWeight = addThis.weight;
35
36 // copy over the first
37 AmgVector(5) finalParameters(firstParameters);
38 double finalWeight = firstWeight;
39 Trk::MultiComponentStateCombiner::combineParametersWithWeight(
40 finalParameters, finalWeight, secondParameters, secondWeight);
41
42 const AmgSymMatrix(5)* firstMeasuredCov = firstTrackParameters->covariance();
43 const AmgSymMatrix(5)* secondMeasuredCov = secondTrackParameters->covariance();
44 // Check to see if first track parameters are measured or not
45 if (firstMeasuredCov && secondMeasuredCov) {
46 AmgSymMatrix(5) finalMeasuredCov(*firstMeasuredCov);
47 Trk::MultiComponentStateCombiner::combineCovWithWeight(
48 firstParameters, finalMeasuredCov, firstWeight, secondParameters,
49 *secondMeasuredCov, secondWeight);
50 mergeTo.params->updateParameters(finalParameters, finalMeasuredCov);
51 mergeTo.weight = finalWeight;
52 } else {
53 mergeTo.params->updateParameters(finalParameters);
54 mergeTo.weight = finalWeight;
55 }
56}
57
59Trk::MultiComponentState mergeFullDistArray(
61 Trk::MultiComponentState&& statesToMerge,
62 const unsigned int maximumNumberOfComponents) {
63
64 const int n = statesToMerge.size();
65 GSFUtils::Component1DArray componentsArray(n);
66 for (int i = 0; i < n; ++i) {
67 const AmgSymMatrix(5)* measuredCov = statesToMerge[i].params->covariance();
68 const AmgVector(5)& parameters = statesToMerge[i].params->parameters();
69 // Fill in infomation
70 const double cov = measuredCov ? (*measuredCov)(Trk::qOverP, Trk::qOverP) : -1.;
71 componentsArray[i].mean = parameters[Trk::qOverP];
72 componentsArray[i].cov = cov;
73 componentsArray[i].invCov = cov > 0 ? 1. / cov : 1e10;
74 componentsArray[i].weight = statesToMerge[i].weight;
75 }
76
77 // Gather the merges
78 const GSFUtils::MergeArray merges = findMerges(std::move(componentsArray), maximumNumberOfComponents);
79
80 // Do the full 5D calculations of the merge
81 const int numMerges = merges.size();
82 for (int i = 0; i < numMerges; ++i) {
83 const int mini = merges[i].To;
84 const int minj = merges[i].From;
85 combineWithWeight(statesToMerge[mini], statesToMerge[minj]);
86 statesToMerge[minj].params.reset();
87 statesToMerge[minj].weight = 0.;
88 }
89 // Assemble the final result
90 for (auto& state : statesToMerge) {
91 // Avoid merge ones
92 if (!state.params) {
93 continue;
94 }
95 cache.multiComponentState.push_back({std::move(state.params), state.weight});
96 cache.validWeightSum += state.weight;
97 }
98 Trk::MultiComponentState mergedState =
100 // Clear the state vector
101 return mergedState;
102}
103
106merge(Trk::MultiComponentState&& statesToMerge,
107 const unsigned int maximumNumberOfComponents) {
108
110 //In case there is nothing to merge
111 if (statesToMerge.size() <= maximumNumberOfComponents) {
112 Trk::MultiComponentStateAssembler::addMultiState(cache, std::move(statesToMerge));
114 }
115
116 // Scan all components for covariance matrices. If one or more component
117 // is missing an error matrix, component reduction is impossible.
118 bool componentWithoutMeasurement = false;
119 Trk::MultiComponentState::const_iterator component = statesToMerge.cbegin();
120 for (; component != statesToMerge.cend(); ++component) {
121 if (!component->params->covariance()) {
122 componentWithoutMeasurement = true;
123 break;
124 }
125 }
126 if (componentWithoutMeasurement) {
127 // Sort to select the one with the largest weight
128 std::sort(
129 statesToMerge.begin(), statesToMerge.end(),
130 [](const Trk::ComponentParameters& x,
131 const Trk::ComponentParameters& y) { return x.weight > y.weight; });
132 //And return it
133 Trk::ComponentParameters dummyCompParams = {std::move(statesToMerge.begin()->params), 1.};
134 Trk::MultiComponentState returnMultiState;
135 returnMultiState.push_back(std::move(dummyCompParams));
136 return returnMultiState;
137 }
138 //Do the full merging of states
139 return mergeFullDistArray(cache, std::move(statesToMerge),maximumNumberOfComponents);
140}
141} // end anonymous namespace
142
143
144// The following does heave use of Eigen
145// for covariance. Avoid out-of-line calls
146// to Eigen
150 const Trk::MultiComponentState& forwardsMultiState,
151 const Trk::MultiComponentState& smootherMultiState,
152 unsigned int maximumNumberOfComponents)
153{
154 auto combinedMultiState = std::make_unique<Trk::MultiComponentState>();
155
156 // Loop over all components in forwards multi-state
157 for (const auto& forwardsComponent : forwardsMultiState) {
158 const AmgSymMatrix(5)* forwardCov = forwardsComponent.params->covariance();
159 // Loop over all components in the smoother multi-state
160 for (const auto& smootherComponent : smootherMultiState) {
161 const AmgSymMatrix(5)* smootherCov = smootherComponent.params->covariance();
162
163 // If not covariances return here.
164 if (!smootherCov && !forwardCov) {
165 return {};
166 }
167 //No forward only smoothed.
168 if (!forwardCov) {
169 Trk::ComponentParameters smootherComponentOnly = {
170 smootherComponent.params->uniqueClone(), smootherComponent.weight};
171 combinedMultiState->push_back(std::move(smootherComponentOnly));
172 continue;
173 }
174 //No smoothed only forward.
175 if (!smootherCov) {
176 Trk::ComponentParameters forwardComponentOnly = {
177 forwardsComponent.params->uniqueClone(), forwardsComponent.weight};
178 combinedMultiState->push_back(std::move(forwardComponentOnly));
179 continue;
180 }
181 /* A comment on te Algebra.
182 * For the covariances P and states X
183 * the most pedagogical presentation found
184 * in textbooks is :
185 *
186 * P_k^s = \left( \left( P_{k|k}^f \right)^{-1}
187 * + \left( P_{k|k-1}^b \right)^{-1}
188 *
189 * X_k^s = P_k^s \left( \left( P_{k|k}^f \right)^{-1} X_{k|k}^f
190 * + \left( P_{k|k-1}^b \right)^{-1} X_{k|k-1}^b \right)
191 *
192 * The f, b, s notation refers to forward, backward, and smoothed
193 * values. The k|k-1 subscript for the backward estimates refers to the
194 * update sequence.
195 *
196 * These can be written (via inversion lemmas) as
197 * P_k^s \left( P_{k|k}^f \right)^{-1} =
198 * I - P_k^s \left( P_{k|k-1}^b \right)^{-1}
199 *
200 * X_k^s = X_{k|k}^f + P_k^s \left( P_{k|k-1}^b \right)^{-1}
201 * \left( X_{k|k-1}^b - X_{k|k}^f \right)
202 *
203 * We use K = P_k^s \left( P_{k|k-1}^b \right)^{-1} below
204 */
205 const AmgVector(5) smootherParams = smootherComponent.params->parameters();
206 const AmgVector(5) forwardParams = forwardsComponent.params->parameters();
207 const AmgSymMatrix(5) summedCovariance = *forwardCov + *smootherCov;
208 const AmgSymMatrix(5) invertedSummedCovariance = summedCovariance.inverse();
209 const AmgSymMatrix(5) K = *forwardCov * invertedSummedCovariance;
210 const AmgVector(5) newParameters = forwardParams + K * (smootherParams - forwardParams);
211 AmgSymMatrix(5) covarianceOfNewParameters = K * (*smootherCov);
212
213 std::unique_ptr<Trk::TrackParameters> combinedTrackParameters =
214 (forwardsComponent.params)
215 ->associatedSurface()
216 .createUniqueTrackParameters(
217 newParameters[Trk::loc1],
218 newParameters[Trk::loc2],
219 newParameters[Trk::phi],
220 newParameters[Trk::theta],
221 newParameters[Trk::qOverP],
222 std::move(covarianceOfNewParameters));
223
224 // Determine the scaling factor for the new weighting.
225 // We need to include the prob for suh difference
226 // via the many-dimensional gaussian pdf.
227 const AmgVector(5) parametersDiff = forwardParams - smootherParams;
228 double const exponent = parametersDiff.transpose() *
229 invertedSummedCovariance * parametersDiff;
230 double const weightScalingFactor = exp(-0.5 * exponent);
231 double const combinedWeight = smootherComponent.weight *
232 forwardsComponent.weight *
233 weightScalingFactor;
234 Trk::ComponentParameters combinedComponent = {
235 std::move(combinedTrackParameters), combinedWeight};
236 combinedMultiState->push_back(std::move(combinedComponent));
237 }
238 } //end of double loop
239
240 // Component reduction on the combined state
241 Trk::MultiComponentState mergedState =
242 merge(std::move(*combinedMultiState), maximumNumberOfComponents);
243 // Before return the weights of the states need to be renormalised to one.
245
246 return mergedState;
247}
#define AmgSymMatrix(dim)
#define AmgVector(rows)
Utilities to facilitate the calculation of the KL divergence/distance between components of the mixtu...
if(febId1==febId2)
#define y
#define x
#define ATH_FLATTEN
MergeArray findMerges(Component1DArray &&componentsIn, const int reducedSize)
Return which components need to be merged.
AlignedDynArray< Component1D, GSFConstants::alignment > Component1DArray
std::vector< Merge > MergeArray
void addMultiState(MultiComponentStateAssembler::Cache &cache, Trk::MultiComponentState &&multiComponentState)
Method to add a new Trk::MultiComponentState to the cached Trk::MultiComponentState object under cons...
MultiComponentState assembledState(MultiComponentStateAssembler::Cache &&cache)
Method to return the cached state object - it performs a reweighting before returning the object base...
void renormaliseState(MultiComponentState &, double norm=1)
Performing renormalisation of total state weighting to one.
Trk::MultiComponentState combine(const Trk::MultiComponentState &forwardsMultiState, const Trk::MultiComponentState &smootherMultiState, unsigned int maximumNumberOfComponents)
Helper to combine forward with smoother MultiComponentStates.
Ensure that the ATLAS eigen extensions are properly loaded.
std::vector< ComponentParameters > MultiComponentState
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ loc2
generic first and second local coordinate
Definition ParamDefs.h:35
@ phi
Definition ParamDefs.h:75
@ loc1
Definition ParamDefs.h:34
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition merge.py:1
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
std::unique_ptr< Trk::TrackParameters > params