CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiGaussianStateTransform.cc
Go to the documentation of this file.
7 
10 {
11  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
12  return multiState(extra->outerStateLocalParameters(),
13  extra->outerStateCovariances(),
14  extra->outerStateWeights());
15 }
16 
19 {
20  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
21  return multiState(extra->innerStateLocalParameters(),
22  extra->innerStateCovariances(),
23  extra->innerStateWeights());
24 }
25 
28  unsigned int index)
29 {
30  if ( index>=N )
31  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
32 
33  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
34  return multiState1D(extra->outerStateLocalParameters(),
35  extra->outerStateCovariances(),
36  extra->outerStateWeights(),
37  index);
38 }
39 
42  unsigned int index)
43 {
44  if ( index>=N )
45  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
46 
47  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
48  return multiState1D(extra->innerStateLocalParameters(),
49  extra->innerStateCovariances(),
50  extra->innerStateWeights(),
51  index);
52 }
53 
56  const std::vector<MultiGaussianState<N>::Matrix>& covariances,
57  const std::vector<double>& weights)
58 {
59  unsigned int nc = parameters.size();
61  components.reserve(nc);
62  for ( unsigned int i=0; i<nc; ++i ) {
64  sgs(new MultiGaussianState<N>::SingleState(parameters[i],covariances[i],weights[i]));
65  components.push_back(sgs);
66  }
68 }
69 
72  const std::vector<MultiGaussianState<N>::Matrix>& covariances,
73  const std::vector<double>& weights, unsigned int index)
74 {
75  unsigned int nc = parameters.size();
77  components.reserve(nc);
78  for ( unsigned int i=0; i<nc; ++i ) {
79  components.push_back(SingleGaussianState1D(parameters[i](index),
80  covariances[i](index,index),
81  weights[i]));
82  }
83  return MultiGaussianState1D(components);
84 }
85 
88 {
89  std::vector<TrajectoryStateOnSurface> tsosComponents(tsos.components());
91  components.reserve(tsosComponents.size());
92  for ( std::vector<TrajectoryStateOnSurface>::const_iterator i=tsosComponents.begin();
93  i!=tsosComponents.end(); ++i ) {
95  sgs(new MultiGaussianState<5>::SingleState(i->localParameters().vector(),
96  i->localError().matrix(),
97  i->weight()));
98  components.push_back(sgs);
99  }
101 }
102 
105  unsigned int index)
106 {
107  if ( index>=N )
108  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
109  std::vector<TrajectoryStateOnSurface> tsosComponents(tsos.components());
111  components.reserve(tsosComponents.size());
112  for ( std::vector<TrajectoryStateOnSurface>::const_iterator i=tsosComponents.begin();
113  i!=tsosComponents.end(); ++i ) {
114  components.push_back(SingleGaussianState1D(i->localParameters().vector()(index),
115  i->localError().matrix()(index,index),
116  i->weight()));
117  }
118  return MultiGaussianState1D(components);
119 }
120 
123  const TrajectoryStateOnSurface refTsos)
124 {
125  const LocalTrajectoryParameters& refPars(refTsos.localParameters());
126  double pzSign = refPars.pzSign();
127  bool charged = refPars.charge()!=0;
128  LocalTrajectoryParameters pars(singleState.mean(),pzSign,charged);
129  LocalTrajectoryError errs(singleState.covariance());
130  // return state (doesn't use weight of the single state)
131  return TrajectoryStateOnSurface(pars,errs,refTsos.surface(),refTsos.magneticField());
132 }
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
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)
std::vector< SingleStatePtr > SingleStateContainer
const Matrix & covariance() const
covariance matrix
const SurfaceType & surface() const
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)
const GsfTrackExtraRef & gsfExtra() const
reference to &quot;extra&quot; object
Definition: GsfTrack.h:32
#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.
std::vector< TrajectoryStateOnSurface > components() const
boost::shared_ptr< SingleState > SingleStatePtr