ATLAS Offline Software
MultiComponentStateCombiner.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
17 //
18 #include "CxxUtils/phihelper.h"
19 #include "CxxUtils/inline_hints.h"
21 #include "TrkSurfaces/Surface.h"
22 
23 namespace {
24 
26 Trk::MultiComponentState mergeFullDistArray(
28  Trk::MultiComponentState&& statesToMerge,
29  const unsigned int maximumNumberOfComponents) {
30  GSFUtils::Component1DArray componentsArray;
31  const int32_t n = statesToMerge.size();
32  componentsArray.numComponents = n;
33  for (int32_t i = 0; i < n; ++i) {
34  const AmgSymMatrix(5)* measuredCov = statesToMerge[i].params->covariance();
35  const AmgVector(5)& parameters = statesToMerge[i].params->parameters();
36  // Fill in infomation
37  const double cov =
38  measuredCov ? (*measuredCov)(Trk::qOverP, Trk::qOverP) : -1.;
39  componentsArray.components[i].mean = parameters[Trk::qOverP];
40  componentsArray.components[i].cov = cov;
41  componentsArray.components[i].invCov = cov > 0 ? 1. / cov : 1e10;
42  componentsArray.components[i].weight = statesToMerge[i].weight;
43  }
44 
45  // Gather the merges
46  const GSFUtils::MergeArray KL =
47  findMerges(componentsArray, maximumNumberOfComponents);
48 
49  // Do the full 5D calculations of the merge
50  const int32_t numMerges = KL.numMerges;
51  for (int32_t i = 0; i < numMerges; ++i) {
52  const int8_t mini = KL.merges[i].To;
53  const int8_t minj = KL.merges[i].From;
55  statesToMerge[minj]);
56  statesToMerge[minj].params.reset();
57  statesToMerge[minj].weight = 0.;
58  }
59  // Assemble the final result
60  for (auto& state : statesToMerge) {
61  // Avoid merge ones
62  if (!state.params) {
63  continue;
64  }
65  cache.multiComponentState.push_back({std::move(state.params),
66  state.weight});
67  cache.validWeightSum += state.weight;
68  }
69  Trk::MultiComponentState mergedState =
71  // Clear the state vector
72  return mergedState;
73 }
74 
76  Trk::MultiComponentState&& statesToMerge,
77  const unsigned int maximumNumberOfComponents) {
78  // Assembler Cache
80  if (statesToMerge.size() <= maximumNumberOfComponents) {
82  std::move(statesToMerge));
84  }
85 
86  // Scan all components for covariance matrices. If one or more component
87  // is missing an error matrix, component reduction is impossible.
88  bool componentWithoutMeasurement = false;
89  Trk::MultiComponentState::const_iterator component = statesToMerge.cbegin();
90  for (; component != statesToMerge.cend(); ++component) {
91  if (!component->params->covariance()) {
92  componentWithoutMeasurement = true;
93  break;
94  }
95  }
96  if (componentWithoutMeasurement) {
97  // Sort to select the one with the largest weight
98  std::sort(statesToMerge.begin(), statesToMerge.end(),
100  return x.weight > y.weight;
101  });
102 
103  Trk::ComponentParameters dummyCompParams = {std::move(statesToMerge.begin()->params), 1.};
104  Trk::MultiComponentState returnMultiState;
105  returnMultiState.push_back(std::move(dummyCompParams));
106  return returnMultiState;
107  }
108 
109  return mergeFullDistArray(cache, std::move(statesToMerge),
110  maximumNumberOfComponents);
111 }
112 // Actual implementation method for combining
113 // a multi component state.
114 Trk::ComponentParameters combineToSingleImpl(
115  const Trk::MultiComponentState& uncombinedState, const bool useMode) {
116  if (uncombinedState.empty()) {
117  return {};
118  }
119 
120  const Trk::TrackParameters* firstParameters =
121  uncombinedState.front().params.get();
122 
123  // Check to see if first track parameters are measured or not
124  const AmgSymMatrix(5)* firstMeasuredCov = firstParameters->covariance();
125 
126  if (uncombinedState.size() == 1) {
127  return {uncombinedState.front().params->uniqueClone(),
128  uncombinedState.front().weight};
129  }
130 
131  double sumW(0.);
132  const int dimension = (uncombinedState.front()).params->parameters().rows();
133 
134  AmgVector(5) mean;
135  mean.setZero();
136  AmgSymMatrix(5) covariance;
137  AmgSymMatrix(5) covariancePart1;
138  covariancePart1.setZero();
139  AmgSymMatrix(5) covariancePart2;
140  covariancePart2.setZero();
141 
142  Trk::MultiComponentState::const_iterator component = uncombinedState.begin();
143  double totalWeight(0.);
144  for (; component != uncombinedState.end(); ++component) {
145  double weight = (*component).weight;
146  totalWeight += weight;
147  }
148 
149  component = uncombinedState.begin();
150 
151  for (; component != uncombinedState.end(); ++component) {
152 
153  const Trk::TrackParameters* trackParameters = (*component).params.get();
154  double weight = (*component).weight;
155 
156  AmgVector(5) parameters = trackParameters->parameters();
157 
158  // Ensure that we don't have any problems with the cyclical nature of phi
159  // Use first state as reference poin
160  // Ensure that we don't have any problems with the cyclical nature of phi
161  // Use first state as reference poin
162  double deltaPhi =
163  (uncombinedState.begin())->params->parameters()[2] - parameters[2];
164  if (deltaPhi > M_PI) {
165  parameters[2] += 2 * M_PI;
166  } else if (deltaPhi < -M_PI) {
167  parameters[2] -= 2 * M_PI;
168  }
169 
170  sumW += weight;
171  mean += weight * parameters;
172 
173  // Extract local error matrix: Must make sure track parameters are measured,
174  // ie have an associated error matrix.
175  const AmgSymMatrix(5)* measuredCov = trackParameters->covariance();
176 
177  // Calculate the combined covariance matrix
178  // \sigma = \Sum_{m=1}^{M} w_{m}(\sigma_m + (\mu_m-\mu)(\mu_m-\mu)^{T})
179  if (measuredCov) {
180  // Changed from errorMatrixInMeasurementFrame
181  covariancePart1 += weight * (*measuredCov);
182  // Loop over all remaining components to find the second part of the
183  // covariance
184  Trk::MultiComponentState::const_iterator remainingComponentIterator =
185  component;
186 
187  for (; remainingComponentIterator != uncombinedState.end();
188  ++remainingComponentIterator) {
189 
190  if (remainingComponentIterator == component) {
191  continue;
192  }
193 
194  AmgVector(5) parameterDifference =
195  parameters - ((*remainingComponentIterator).params)->parameters();
196 
197  double remainingComponentIteratorWeight =
198  (*remainingComponentIterator).weight;
199 
200  covariancePart2 += weight * remainingComponentIteratorWeight *
201  parameterDifference *
202  parameterDifference.transpose();
203 
204  } // end loop over remaining components
205  } // end clause if errors are involved
206 
207  if (weight / totalWeight > 1.0) {
208  break;
209  }
210  } // end loop over all components
211 
212  mean /= sumW;
213 
214  // Ensure that phi is between -pi and pi
215  //
216  mean[2] = CxxUtils::wrapToPi(mean[2]);
217 
218  covariance = covariancePart1 / sumW + covariancePart2 / (sumW * sumW);
219 
220  if (useMode && dimension == 5) {
221 
222  // Calculate the mode of the q/p distribution
223  std::array<double, 10> modes =
225 
226  // Replace mean with mode if qOverP mode is not 0
227  if (modes[4] != 0) {
228  mean[0] = modes[0];
229  mean[1] = modes[1];
230  mean[2] = modes[2];
231  mean[3] = modes[3];
232  mean[4] = modes[4];
233 
234  if (modes[5 + 0] > 0) {
235  double currentErr = sqrt((covariance)(0, 0));
236  currentErr = modes[5 + 0] / currentErr;
237  (covariance)(0, 0) = modes[5 + 0] * modes[5 + 0];
238  covariance.fillSymmetric(1, 0, (covariance)(1, 0) * currentErr);
239  covariance.fillSymmetric(2, 0, (covariance)(2, 0) * currentErr);
240  covariance.fillSymmetric(3, 0, (covariance)(3, 0) * currentErr);
241  covariance.fillSymmetric(4, 0, (covariance)(4, 0) * currentErr);
242  }
243  if (modes[5 + 1] > 0) {
244  double currentErr = sqrt((covariance)(1, 1));
245  currentErr = modes[5 + 1] / currentErr;
246  covariance.fillSymmetric(1, 0, (covariance)(1, 0) * currentErr);
247  (covariance)(1, 1) = modes[5 + 1] * modes[5 + 1];
248  covariance.fillSymmetric(2, 1, (covariance)(2, 1) * currentErr);
249  covariance.fillSymmetric(3, 1, (covariance)(3, 1) * currentErr);
250  covariance.fillSymmetric(4, 1, (covariance)(4, 1) * currentErr);
251  }
252  if (modes[5 + 2] > 0) {
253  double currentErr = sqrt((covariance)(2, 2));
254  currentErr = modes[5 + 2] / currentErr;
255  covariance.fillSymmetric(2, 0, (covariance)(2, 0) * currentErr);
256  covariance.fillSymmetric(2, 1, (covariance)(2, 1) * currentErr);
257  (covariance)(2, 2) = modes[5 + 2] * modes[5 + 2];
258  covariance.fillSymmetric(3, 2, (covariance)(3, 2) * currentErr);
259  covariance.fillSymmetric(4, 2, (covariance)(4, 2) * currentErr);
260  }
261  if (modes[5 + 3] > 0) {
262  double currentErr = sqrt((covariance)(3, 3));
263  currentErr = modes[5 + 3] / currentErr;
264  covariance.fillSymmetric(3, 0, (covariance)(3, 0) * currentErr);
265  covariance.fillSymmetric(3, 1, (covariance)(3, 1) * currentErr);
266  covariance.fillSymmetric(3, 2, (covariance)(3, 2) * currentErr);
267  (covariance)(3, 3) = modes[5 + 3] * modes[5 + 3];
268  covariance.fillSymmetric(4, 3, (covariance)(4, 3) * currentErr);
269  }
270  if (modes[5 + 4] > 0) {
271  double currentErr = sqrt((covariance)(4, 4));
272  currentErr = modes[5 + 4] / currentErr;
273  covariance.fillSymmetric(4, 0, (covariance)(4, 0) * currentErr);
274  covariance.fillSymmetric(4, 1, (covariance)(4, 1) * currentErr);
275  covariance.fillSymmetric(4, 2, (covariance)(4, 2) * currentErr);
276  covariance.fillSymmetric(4, 3, (covariance)(4, 3) * currentErr);
277  (covariance)(4, 4) = modes[5 + 4] * modes[5 + 4];
278  }
279 
280  } // modes[4]!=0
281  } // useMode && dimensions==5
282 
283  std::unique_ptr<Trk::TrackParameters> combinedTrackParameters = nullptr;
284  double loc1 = mean[Trk::loc1];
285  double loc2 = mean[Trk::loc2];
286  double phi = mean[Trk::phi];
287  double theta = mean[Trk::theta];
288  double qoverp = mean[Trk::qOverP];
289  if (firstMeasuredCov) {
290  combinedTrackParameters =
292  loc1, loc2, phi, theta, qoverp, std::move(covariance));
293  } else {
294  combinedTrackParameters =
296  loc1, loc2, phi, theta, qoverp, std::nullopt);
297  }
298 
299  return {std::move(combinedTrackParameters), totalWeight};
300 }
301 
302 } // end anonymous namespace
303 
304 std::unique_ptr<Trk::TrackParameters>
306 const Trk::MultiComponentState& uncombinedState, const bool useMode) {
307  Trk::ComponentParameters combinedComponent =
308  combineToSingleImpl(uncombinedState, useMode);
309  return std::move(combinedComponent.params);
310 }
311 
312 void
314 Trk::ComponentParameters& mergeTo,
315 const Trk::ComponentParameters& addThis)
316 {
317  const Trk::TrackParameters* firstTrackParameters = mergeTo.params.get();
318  const AmgVector(5)& firstParameters = firstTrackParameters->parameters();
319  double firstWeight = mergeTo.weight;
320 
321  const Trk::TrackParameters* secondTrackParameters = addThis.params.get();
322  const AmgVector(5)& secondParameters = secondTrackParameters->parameters();
323  double secondWeight = addThis.weight;
324 
325  // copy over the first
326  AmgVector(5) finalParameters(firstParameters);
327  double finalWeight = firstWeight;
329  finalParameters, finalWeight, secondParameters, secondWeight);
330 
331  const AmgSymMatrix(5)* firstMeasuredCov = firstTrackParameters->covariance();
332  const AmgSymMatrix(5)* secondMeasuredCov = secondTrackParameters->covariance();
333  // Check to see if first track parameters are measured or not
334  if (firstMeasuredCov && secondMeasuredCov) {
335  AmgSymMatrix(5) finalMeasuredCov(*firstMeasuredCov);
336  combineCovWithWeight(firstParameters, finalMeasuredCov, firstWeight,
337  secondParameters, *secondMeasuredCov, secondWeight);
338  mergeTo.params->updateParameters(finalParameters, finalMeasuredCov);
339  mergeTo.weight = finalWeight;
340  } else {
341  mergeTo.params->updateParameters(finalParameters);
342  mergeTo.weight = finalWeight;
343  }
344 }
353 // The following does heave use of Eigen
354 // for covariance. Avoid out-of-line calls
355 // to Eigen
357 void
359  AmgVector(5) & firstParameters,
360  double& firstWeight,
361  const AmgVector(5) & secondParameters,
362  const double secondWeight)
363 {
364  double totalWeight = firstWeight + secondWeight;
365  double invTotalWeight = 1.0/totalWeight;
366  double deltaPhi = firstParameters[2] - secondParameters[2];
367  if (deltaPhi > M_PI) {
368  firstParameters[2] -= 2 * M_PI;
369  } else if (deltaPhi < -M_PI) {
370  firstParameters[2] += 2 * M_PI;
371  }
372  firstParameters =
373  (firstWeight * firstParameters + secondWeight * secondParameters) *
374  invTotalWeight;
375  // Ensure that phi is between -pi and pi
376  firstParameters[2] = CxxUtils::wrapToPi(firstParameters[2]);
377  firstWeight = totalWeight;
378 }
379 
380 // The following does heave use of Eigen
381 // for covariance. Avoid out-of-line calls
382 // to Eigen
384 void
386  const AmgVector(5) & firstParameters,
387  AmgSymMatrix(5) & firstMeasuredCov,
388  const double firstWeight,
389  const AmgVector(5) & secondParameters,
390  const AmgSymMatrix(5) & secondMeasuredCov,
391  const double secondWeight)
392 {
393  double invTotalWeight = 1.0/(firstWeight + secondWeight);
394  AmgVector(5) parameterDifference = firstParameters - secondParameters;
395  parameterDifference[2] = CxxUtils::wrapToPi(parameterDifference[2]);
396  parameterDifference *= invTotalWeight;
397  firstMeasuredCov =
398  (firstWeight * firstMeasuredCov + secondWeight * secondMeasuredCov) *
399  invTotalWeight;
400  firstMeasuredCov += firstWeight * secondWeight * parameterDifference *
401  parameterDifference.transpose();
402 }
403 
404 // The following does heave use of Eigen
405 // for covariance. Avoid out-of-line calls
406 // to Eigen
410  const Trk::MultiComponentState& forwardsMultiState,
411  const Trk::MultiComponentState& smootherMultiState,
412  unsigned int maximumNumberOfComponents)
413 {
414 
415  std::unique_ptr<Trk::MultiComponentState> combinedMultiState =
416  std::make_unique<Trk::MultiComponentState>();
417 
418  // Loop over all components in forwards multi-state
419  for (const auto& forwardsComponent : forwardsMultiState) {
420  // Need to check that all components have associated weight matricies
421  const AmgSymMatrix(5)* forwardMeasuredCov =
422  forwardsComponent.params->covariance();
423  // Loop over all components in the smoother multi-state
424  for (const auto& smootherComponent : smootherMultiState) {
425  // Need to check that all components have associated weight matricies
426  const AmgSymMatrix(5)* smootherMeasuredCov =
427  smootherComponent.params->covariance();
428  if (!smootherMeasuredCov && !forwardMeasuredCov) {
429  return {};
430  }
431 
432  if (!forwardMeasuredCov) {
433  Trk::ComponentParameters smootherComponentOnly = {
434  smootherComponent.params->uniqueClone(), smootherComponent.weight};
435  combinedMultiState->push_back(std::move(smootherComponentOnly));
436  continue;
437  }
438 
439  if (!smootherMeasuredCov) {
440  Trk::ComponentParameters forwardComponentOnly = {
441  forwardsComponent.params->uniqueClone(), forwardsComponent.weight};
442  combinedMultiState->push_back(std::move(forwardComponentOnly));
443  continue;
444  }
445 
446  const AmgSymMatrix(5) summedCovariance =
447  *forwardMeasuredCov + *smootherMeasuredCov;
448  const AmgSymMatrix(5) K =
449  *forwardMeasuredCov * summedCovariance.inverse();
450  const AmgVector(5) newParameters =
451  forwardsComponent.params->parameters() +
452  K * (smootherComponent.params->parameters() -
453  forwardsComponent.params->parameters());
454  const AmgVector(5) parametersDiff =
455  forwardsComponent.params->parameters() -
456  smootherComponent.params->parameters();
457 
458  AmgSymMatrix(5) covarianceOfNewParameters =
459  AmgSymMatrix(5)(K * *smootherMeasuredCov);
460 
461  std::unique_ptr<Trk::TrackParameters> combinedTrackParameters =
462  (forwardsComponent.params)
463  ->associatedSurface()
464  .createUniqueTrackParameters(newParameters[Trk::loc1],
465  newParameters[Trk::loc2],
466  newParameters[Trk::phi],
467  newParameters[Trk::theta],
468  newParameters[Trk::qOverP],
469  std::move(covarianceOfNewParameters));
470  const AmgSymMatrix(5) invertedSummedCovariance =
471  summedCovariance.inverse();
472  // Determine the scaling factor for the new weighting. Determined from the
473  // PDF of the many-dimensional gaussian
474  double exponent =
475  parametersDiff.transpose() * invertedSummedCovariance * parametersDiff;
476  double weightScalingFactor = exp(-0.5 * exponent);
477  double combinedWeight = smootherComponent.weight *
478  forwardsComponent.weight * weightScalingFactor;
479  Trk::ComponentParameters combinedComponent = {
480  std::move(combinedTrackParameters), combinedWeight};
481  combinedMultiState->push_back(std::move(combinedComponent));
482  }
483  }
484  // Component reduction on the combined state
485  Trk::MultiComponentState mergedState =
486  merge(std::move(*combinedMultiState), maximumNumberOfComponents);
487  // Before return the weights of the states need to be renormalised to one.
489 
490  return mergedState;
491 }
MultiComponentStateCombiner.h
mean
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
Definition: dependence.cxx:254
inline_hints.h
Trk::ComponentParameters::weight
double weight
Definition: ComponentParameters.h:24
GSFUtils::MergeArray
struct representing an array or the merges.
Definition: KLGaussianMixtureReduction.h:54
TrackParameters.h
Trk::MultiComponentStateAssembler::assembledState
MultiComponentState assembledState(MultiComponentStateAssembler::Cache &&cache)
Method to return the cached state object - it performs a reweighting before returning the object base...
Definition: MultiComponentStateAssembler.cxx:120
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:64
Surface.h
Trk::ParametersBase::associatedSurface
virtual const Surface & associatedSurface() const override=0
Access to the Surface associated to the Parameters.
Trk::ParametersBase::uniqueClone
std::unique_ptr< ParametersBase< DIM, T > > uniqueClone() const
clone method for polymorphic deep copy returning unique_ptr; it is not overriden, but uses the existi...
Definition: ParametersBase.h:97
xAOD::deltaPhi
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap setEtaBin setIsTgcFailure setDeltaPt deltaPhi
Definition: L2StandAloneMuon_v1.cxx:160
CxxUtils::wrapToPi
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition: phihelper.h:24
theta
Scalar theta() const
theta method
Definition: AmgMatrixBasePlugin.h:71
PlotCalibFromCool.begin
begin
Definition: PlotCalibFromCool.py:94
plotBeamSpotVxVal.cov
cov
Definition: plotBeamSpotVxVal.py:201
M_PI
#define M_PI
Definition: ActiveFraction.h:11
Trk::loc2
@ loc2
generic first and second local coordinate
Definition: ParamDefs.h:41
Trk::MultiComponentStateCombiner::combineToSingle
std::unique_ptr< Trk::TrackParameters > combineToSingle(const MultiComponentState &, const bool useMode=false)
@bried Calculate combined state of many components
Definition: MultiComponentStateCombiner.cxx:305
Trk::MultiComponentStateAssembler::Cache
Definition: MultiComponentStateAssembler.h:35
drawFromPickle.exp
exp
Definition: drawFromPickle.py:36
x
#define x
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:52
Trk::MultiComponentStateModeCalculator::calculateMode
std::array< double, 10 > calculateMode(const MultiComponentState &)
Method to calculate mode with MultiComponentState state as input.
Definition: MultiComponentStateModeCalculator.cxx:386
mergePhysValFiles.end
end
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:93
Trk::MultiComponentStateCombiner::combineCovWithWeight
void combineCovWithWeight(const AmgVector(5) &firstParameters, AmgSymMatrix(5) &firstMeasuredCov, const double firstWeight, const AmgVector(5) &secondParameters, const AmgSymMatrix(5) &secondMeasuredCov, const double secondWeight)
Update cov matrix.
Definition: MultiComponentStateCombiner.cxx:385
dqt_zlumi_pandas.weight
int weight
Definition: dqt_zlumi_pandas.py:200
GSFUtils::MergeArray::numMerges
int32_t numMerges
Definition: KLGaussianMixtureReduction.h:61
MultiComponentStateAssembler.h
lumiFormat.i
int i
Definition: lumiFormat.py:92
Trk::Surface::createUniqueTrackParameters
virtual ChargedTrackParametersUniquePtr createUniqueTrackParameters(double l1, double l2, double phi, double theat, double qop, std::optional< AmgSymMatrix(5)> cov=std::nullopt) const =0
Use the Surface as a ParametersBase constructor, from local parameters - charged.
ATH_FLATTEN
#define ATH_FLATTEN
Definition: inline_hints.h:52
beamspotman.n
n
Definition: beamspotman.py:731
Trk::theta
@ theta
Definition: ParamDefs.h:72
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
Trk::MultiComponentState
std::vector< ComponentParameters > MultiComponentState
Definition: ComponentParameters.h:27
MultiComponentStateModeCalculator.h
Trk::ParametersBase
Definition: ParametersBase.h:55
Trk::MultiComponentStateCombiner::combineWithWeight
void combineWithWeight(Trk::ComponentParameters &mergeTo, const Trk::ComponentParameters &addThis)
Combined/merge a component to another one.
Definition: MultiComponentStateCombiner.cxx:313
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Trk::MultiComponentStateCombiner::combineParametersWithWeight
void combineParametersWithWeight(AmgVector(5) &firstParameters, double &firstWeight, const AmgVector(5) &secondParameters, const double secondWeight)
Update parameters.
Definition: MultiComponentStateCombiner.cxx:358
Trk::MultiComponentStateHelpers::renormaliseState
void renormaliseState(MultiComponentState &, double norm=1)
Performing renormalisation of total state weighting to one.
Definition: ComponentParameters.cxx:80
phihelper.h
Helper for azimuthal angle calculations.
GSFUtils::Component1DArray::components
std::array< Component1D, GSFConstants::maxComponentsAfterConvolution > components
Definition: KLGaussianMixtureReduction.h:45
Trk::ComponentParameters
Definition: ComponentParameters.h:22
KLGaussianMixtureReduction.h
Utilities to facilitate the calculation of the KL divergence/distance between components of the mixtu...
y
#define y
Trk::MultiComponentStateAssembler::Cache::multiComponentState
Trk::MultiComponentState multiComponentState
Definition: MultiComponentStateAssembler.h:39
GSFUtils::Component1DArray::numComponents
int32_t numComponents
Definition: KLGaussianMixtureReduction.h:46
std::sort
void sort(typename std::reverse_iterator< DataModel_detail::iterator< DVL > > beg, typename std::reverse_iterator< DataModel_detail::iterator< DVL > > end, const Compare &comp)
Specialization of sort for DataVector/List.
Definition: DVL_algorithms.h:623
Trk::qOverP
@ qOverP
perigee
Definition: ParamDefs.h:73
if
if(febId1==febId2)
Definition: LArRodBlockPhysicsV0.cxx:569
python.utility.LHE.merge
def merge(input_file_pattern, output_file)
Merge many input LHE files into a single output file.
Definition: LHE.py:17
Trk::MultiComponentStateAssembler::addMultiState
void addMultiState(MultiComponentStateAssembler::Cache &cache, Trk::MultiComponentState &&multiComponentState)
Method to add a new Trk::MultiComponentState to the cached Trk::MultiComponentState object under cons...
Definition: MultiComponentStateAssembler.h:79
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
Trk::phi
@ phi
Definition: ParamDefs.h:81
Trk::MultiComponentStateCombiner::combineWithSmoother
Trk::MultiComponentState combineWithSmoother(const Trk::MultiComponentState &forwardsMultiState, const Trk::MultiComponentState &smootherMultiState, unsigned int maximumNumberOfComponents)
Helper to combine forward with smoother MultiComponentStates.
Definition: MultiComponentStateCombiner.cxx:409
PowhegControl_ttFCNC_NLO.params
params
Definition: PowhegControl_ttFCNC_NLO.py:226
Trk::loc1
@ loc1
Definition: ParamDefs.h:40
Trk::MultiComponentStateAssembler::Cache::validWeightSum
double validWeightSum
Definition: MultiComponentStateAssembler.h:40
Trk::ComponentParameters::params
std::unique_ptr< Trk::TrackParameters > params
Definition: ComponentParameters.h:23
GSFUtils::findMerges
MergeArray findMerges(const Component1DArray &componentsIn, const int8_t reducedSize)
Find the order in which the components need to be merged.
Definition: KLGaussianMixtureReduction.cxx:394
GSFUtils::Component1DArray
struct representing an array of 1D component.
Definition: KLGaussianMixtureReduction.h:42
GSFUtils::MergeArray::merges
std::array< merge, GSFConstants::maxComponentsAfterConvolution > merges
Definition: KLGaussianMixtureReduction.h:60