CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MultiGaussianStateTransform.cc
Go to the documentation of this file.
8 
10  const reco::GsfTrack& tk) {
11  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
12  return multiState(extra->outerStateLocalParameters(), extra->outerStateCovariances(), extra->outerStateWeights());
13 }
14 
16  const reco::GsfTrack& tk) {
17  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
18  return multiState(extra->innerStateLocalParameters(), extra->innerStateCovariances(), extra->innerStateWeights());
19 }
20 
22  if (index >= N)
23  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
24 
25  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
26  return multiState1D(
27  extra->outerStateLocalParameters(), extra->outerStateCovariances(), extra->outerStateWeights(), index);
28 }
29 
31  if (index >= N)
32  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
33 
34  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
35  return multiState1D(
36  extra->innerStateLocalParameters(), extra->innerStateCovariances(), extra->innerStateWeights(), index);
37 }
38 
41  const std::vector<MultiGaussianState<N>::Matrix>& covariances,
42  const std::vector<double>& weights) {
43  unsigned int nc = parameters.size();
45  components.reserve(nc);
46  for (unsigned int i = 0; i < nc; ++i) {
48  new MultiGaussianState<N>::SingleState(parameters[i], covariances[i], weights[i]));
49  components.push_back(sgs);
50  }
52 }
53 
56  const std::vector<MultiGaussianState<N>::Matrix>& covariances,
57  const std::vector<double>& weights,
58  unsigned int index) {
59  unsigned int nc = parameters.size();
61  components.reserve(nc);
62  for (unsigned int i = 0; i < nc; ++i) {
63  components.push_back(SingleGaussianState1D(parameters[i](index), covariances[i](index, index), weights[i]));
64  }
65  return MultiGaussianState1D(components);
66 }
67 
69  GetComponents comps(tsos);
70  auto const& tsosComponents = comps();
72  components.reserve(tsosComponents.size());
73  for (auto i = tsosComponents.begin(); i != tsosComponents.end(); ++i) {
75  new MultiGaussianState<5>::SingleState(i->localParameters().vector(), i->localError().matrix(), i->weight()));
76  components.push_back(sgs);
77  }
79 }
80 
82  unsigned int index) {
83  if (index >= N)
84  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
85  GetComponents comps(tsos);
86  auto const& tsosComponents = comps();
88  components.reserve(tsosComponents.size());
89  for (auto i = tsosComponents.begin(); i != tsosComponents.end(); ++i) {
90  components.push_back(SingleGaussianState1D(
91  i->localParameters().vector()(index), i->localError().matrix()(index, index), i->weight()));
92  }
93  return MultiGaussianState1D(components);
94 }
95 
97  const TrajectoryStateOnSurface refTsos) {
98  const LocalTrajectoryParameters& refPars(refTsos.localParameters());
99  double pzSign = refPars.pzSign();
100  bool charged = refPars.charge() != 0;
101  LocalTrajectoryParameters pars(singleState.mean(), pzSign, charged);
102  LocalTrajectoryError errs(singleState.covariance());
103  // return state (doesn't use weight of the single state)
104  return TrajectoryStateOnSurface(pars, errs, refTsos.surface(), refTsos.magneticField());
105 }
MultiGaussianState1D outerMultiState1D(const reco::GsfTrack &tk, unsigned int index)
const LocalTrajectoryParameters & localParameters() const
TrajectoryStateOnSurface tsosFromSingleState(const SingleGaussianState< 5 > &, const TrajectoryStateOnSurface)
const Vector & mean() const
parameter vector
std::vector< SingleGaussianState1D > SingleState1dContainer
const MagneticField * magneticField() const
Mixture of multi-variate gaussian states.
MultiGaussianState< N > innerMultiState(const reco::GsfTrack &tk)
const Matrix & covariance() const
covariance matrix
const SurfaceType & surface() const
std::vector< SingleStatePtr > SingleStateContainer
MultiGaussianState1D multiState1D(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &, unsigned int)
MultiGaussianState< N > outerMultiState(const reco::GsfTrack &tk)
std::shared_ptr< SingleState > SingleStatePtr
const GsfTrackExtraRef & gsfExtra() const
reference to &quot;extra&quot; object
Definition: GsfTrack.h:31
#define N
Definition: blowfish.cc:9
MultiGaussianState< N > multiState(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &)
MultiGaussianState1D innerMultiState1D(const reco::GsfTrack &tk, unsigned int index)
SingleGaussianState< N >::Vector Vector
SingleGaussianState< N >::Matrix Matrix
float pzSign() const
Sign of the z-component of the momentum in the local frame.