ATLAS Offline Software
Loading...
Searching...
No Matches
Trk::TwoStateCombiner Namespace Reference

Functions

Trk::MultiComponentState combine (const Trk::MultiComponentState &forwardsMultiState, const Trk::MultiComponentState &smootherMultiState, unsigned int maximumNumberOfComponents)
 Helper to combine forward with smoother MultiComponentStates.

Function Documentation

◆ combine()

ATH_FLATTEN Trk::MultiComponentState Trk::TwoStateCombiner::combine ( const Trk::MultiComponentState & forwardsMultiState,
const Trk::MultiComponentState & smootherMultiState,
unsigned int maximumNumberOfComponents )

Helper to combine forward with smoother MultiComponentStates.

Definition at line 149 of file TwoStateCombiner.cxx.

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)
void renormaliseState(MultiComponentState &, double norm=1)
Performing renormalisation of total state weighting to one.
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
merge(input_file_pattern, output_file)
Merge many input LHE files into a single output file.
Definition LHE.py:29