CMS 3D CMS Logo

Enumerations | Functions
MultiGaussianStateTransform Namespace Reference

Enumerations

enum  { N = reco::GsfTrackExtra::dimension }
 

Functions

MultiGaussianState< NinnerMultiState (const reco::GsfTrack &tk)
 
MultiGaussianState1D innerMultiState1D (const reco::GsfTrack &tk, unsigned int index)
 
MultiGaussianState< NmultiState (const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &)
 
MultiGaussianState< 5 > multiState (const TrajectoryStateOnSurface)
 
MultiGaussianState1D multiState1D (const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &, unsigned int)
 
MultiGaussianState1D multiState1D (const TrajectoryStateOnSurface, unsigned int)
 
MultiGaussianState< NouterMultiState (const reco::GsfTrack &tk)
 
MultiGaussianState1D outerMultiState1D (const reco::GsfTrack &tk, unsigned int index)
 
TrajectoryStateOnSurface tsosFromSingleState (const SingleGaussianState< 5 > &, const TrajectoryStateOnSurface)
 

Enumeration Type Documentation

anonymous enum
Enumerator

Definition at line 19 of file MultiGaussianStateTransform.h.

Function Documentation

MultiGaussianState< MultiGaussianStateTransform::N > MultiGaussianStateTransform::innerMultiState ( const reco::GsfTrack tk)

Construct a MultiGaussianState from the reco::GsfTrack innermost state (local parameters)

Definition at line 19 of file MultiGaussianStateTransform.cc.

References reco::GsfTrack::gsfExtra(), and multiState().

20 {
21  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
22  return multiState(extra->innerStateLocalParameters(),
23  extra->innerStateCovariances(),
24  extra->innerStateWeights());
25 }
const GsfTrackExtraRef & gsfExtra() const
reference to "extra" object
Definition: GsfTrack.h:32
MultiGaussianState< N > multiState(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &)
MultiGaussianState1D MultiGaussianStateTransform::innerMultiState1D ( const reco::GsfTrack tk,
unsigned int  index 
)

Construct a MultiGaussianState1D for the local parameter corresponding to "index" (0<=index<5) from the reco::GsfTrack innermost state

Definition at line 42 of file MultiGaussianStateTransform.cc.

References Exception, reco::GsfTrack::gsfExtra(), diffTreeTool::index, multiState1D(), and N.

44 {
45  if ( index>=N )
46  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
47 
48  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
49  return multiState1D(extra->innerStateLocalParameters(),
50  extra->innerStateCovariances(),
51  extra->innerStateWeights(),
52  index);
53 }
MultiGaussianState1D multiState1D(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &, unsigned int)
const GsfTrackExtraRef & gsfExtra() const
reference to "extra" object
Definition: GsfTrack.h:32
#define N
Definition: blowfish.cc:9
MultiGaussianState< MultiGaussianStateTransform::N > MultiGaussianStateTransform::multiState ( const std::vector< MultiGaussianState< N >::Vector > &  parameters,
const std::vector< MultiGaussianState< N >::Matrix > &  covariances,
const std::vector< double > &  weights 
)

Construct a MultiGaussianState from the vectors of parameters, covariances and weights

Definition at line 56 of file MultiGaussianStateTransform.cc.

References makeMuonMisalignmentScenario::components, i, and metProducer_cfi::parameters.

Referenced by innerMultiState(), MultiTrajectoryStateMerger::merge(), outerMultiState(), and TagProbeFitter::setWeightVar().

59 {
60  unsigned int nc = parameters.size();
62  components.reserve(nc);
63  for ( unsigned int i=0; i<nc; ++i ) {
65  sgs(new MultiGaussianState<N>::SingleState(parameters[i],covariances[i],weights[i]));
66  components.push_back(sgs);
67  }
69 }
int i
Definition: DBlmapReader.cc:9
Mixture of multi-variate gaussian states.
std::vector< SingleStatePtr > SingleStateContainer
std::shared_ptr< SingleState > SingleStatePtr
MultiGaussianState< 5 > MultiGaussianStateTransform::multiState ( const TrajectoryStateOnSurface  tsos)

Construct a MultiGaussianState from a TrajectoryStateOnSurface (local parameters)

Definition at line 88 of file MultiGaussianStateTransform.cc.

References makeMuonMisalignmentScenario::components, and i.

89 {
90  GetComponents comps(tsos);
91  auto const & tsosComponents = comps();
93  components.reserve(tsosComponents.size());
94  for (auto i=tsosComponents.begin();
95  i!=tsosComponents.end(); ++i ) {
97  sgs(new MultiGaussianState<5>::SingleState(i->localParameters().vector(),
98  i->localError().matrix(),
99  i->weight()));
100  components.push_back(sgs);
101  }
103 }
int i
Definition: DBlmapReader.cc:9
Mixture of multi-variate gaussian states.
std::vector< SingleStatePtr > SingleStateContainer
std::shared_ptr< SingleState > SingleStatePtr
MultiGaussianState1D MultiGaussianStateTransform::multiState1D ( const std::vector< MultiGaussianState< N >::Vector > &  parameters,
const std::vector< MultiGaussianState< N >::Matrix > &  covariances,
const std::vector< double > &  weights,
unsigned int  index 
)

Construct a MultiGaussianState1D from the vectors of parameters, covariances and weights

Definition at line 72 of file MultiGaussianStateTransform.cc.

References makeMuonMisalignmentScenario::components, i, and metProducer_cfi::parameters.

Referenced by MultiTrajectoryStateMode::chargeFromMode(), PFGsfHelper::computeQpMode(), ElectronMomentumCorrector::correct(), GsfTrackProducerBase::fillMode(), innerMultiState1D(), GsfTrackProducerBase::localParametersFromQpMode(), MultiTrajectoryStateMode::momentumFromModeLocal(), MultiTrajectoryStateMode::momentumFromModeP(), MultiTrajectoryStateMode::momentumFromModeQP(), outerMultiState1D(), and MultiTrajectoryStateMode::positionFromModeLocal().

75 {
76  unsigned int nc = parameters.size();
78  components.reserve(nc);
79  for ( unsigned int i=0; i<nc; ++i ) {
80  components.push_back(SingleGaussianState1D(parameters[i](index),
81  covariances[i](index,index),
82  weights[i]));
83  }
84  return MultiGaussianState1D(components);
85 }
int i
Definition: DBlmapReader.cc:9
std::vector< SingleGaussianState1D > SingleState1dContainer
MultiGaussianState1D MultiGaussianStateTransform::multiState1D ( const TrajectoryStateOnSurface  tsos,
unsigned int  index 
)

Construct a MultiGaussianState1D from a TrajectoryStateOnSurface (local parameters)

Definition at line 106 of file MultiGaussianStateTransform.cc.

References makeMuonMisalignmentScenario::components, Exception, i, diffTreeTool::index, and N.

108 {
109  if ( index>=N )
110  throw cms::Exception("LogicError") << "MultiGaussianStateTransform: index out of range";
111  GetComponents comps(tsos);
112  auto const & tsosComponents = comps();
114  components.reserve(tsosComponents.size());
115  for (auto i=tsosComponents.begin();
116  i!=tsosComponents.end(); ++i ) {
117  components.push_back(SingleGaussianState1D(i->localParameters().vector()(index),
118  i->localError().matrix()(index,index),
119  i->weight()));
120  }
121  return MultiGaussianState1D(components);
122 }
int i
Definition: DBlmapReader.cc:9
std::vector< SingleGaussianState1D > SingleState1dContainer
#define N
Definition: blowfish.cc:9
MultiGaussianState< MultiGaussianStateTransform::N > MultiGaussianStateTransform::outerMultiState ( const reco::GsfTrack tk)

Construct a MultiGaussianState from the reco::GsfTrack innermost state (local parameters)

Definition at line 10 of file MultiGaussianStateTransform.cc.

References reco::GsfTrack::gsfExtra(), and multiState().

11 {
12  const reco::GsfTrackExtraRef& extra(tk.gsfExtra());
13  return multiState(extra->outerStateLocalParameters(),
14  extra->outerStateCovariances(),
15  extra->outerStateWeights());
16 }
const GsfTrackExtraRef & gsfExtra() const
reference to "extra" object
Definition: GsfTrack.h:32
MultiGaussianState< N > multiState(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &)
MultiGaussianState1D MultiGaussianStateTransform::outerMultiState1D ( const reco::GsfTrack tk,
unsigned int  index 
)

Construct a MultiGaussianState1D for the local parameter corresponding to "index" (0<=index<5) from the reco::GsfTrack outermost state

Definition at line 28 of file MultiGaussianStateTransform.cc.

References Exception, reco::GsfTrack::gsfExtra(), diffTreeTool::index, multiState1D(), and N.

30 {
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(extra->outerStateLocalParameters(),
36  extra->outerStateCovariances(),
37  extra->outerStateWeights(),
38  index);
39 }
MultiGaussianState1D multiState1D(const std::vector< MultiGaussianState< N >::Vector > &, const std::vector< MultiGaussianState< N >::Matrix > &, const std::vector< double > &, unsigned int)
const GsfTrackExtraRef & gsfExtra() const
reference to "extra" object
Definition: GsfTrack.h:32
#define N
Definition: blowfish.cc:9
TrajectoryStateOnSurface MultiGaussianStateTransform::tsosFromSingleState ( const SingleGaussianState< 5 > &  singleState,
const TrajectoryStateOnSurface  refTsos 
)

Construct a TrajectoryStateOnSurface from a 5D SingleGaussianState (local parameters) and a reference TSOS (surface, charge, ..)

Definition at line 125 of file MultiGaussianStateTransform.cc.

References SingleGaussianState< N >::covariance(), TrajectoryStateOnSurface::localParameters(), TrajectoryStateOnSurface::magneticField(), SingleGaussianState< N >::mean(), LocalTrajectoryParameters::pzSign(), and TrajectoryStateOnSurface::surface().

127 {
128  const LocalTrajectoryParameters& refPars(refTsos.localParameters());
129  double pzSign = refPars.pzSign();
130  bool charged = refPars.charge()!=0;
131  LocalTrajectoryParameters pars(singleState.mean(),pzSign,charged);
132  LocalTrajectoryError errs(singleState.covariance());
133  // return state (doesn't use weight of the single state)
134  return TrajectoryStateOnSurface(pars,errs,refTsos.surface(),refTsos.magneticField());
135 }
const LocalTrajectoryParameters & localParameters() const
const Vector & mean() const
parameter vector
const MagneticField * magneticField() const
const Matrix & covariance() const
covariance matrix
const SurfaceType & surface() const
float pzSign() const
Sign of the z-component of the momentum in the local frame.