CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Protected Member Functions | Protected Attributes
DualKalmanFactory Class Reference
Inheritance diagram for DualKalmanFactory:
TrajectoryFactoryBase

Classes

struct  DualKalmanInput
 

Public Member Functions

virtual DualKalmanFactoryclone () const
 
 DualKalmanFactory (const edm::ParameterSet &config)
 
virtual const
ReferenceTrajectoryCollection 
trajectories (const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const
 Produce the reference trajectories. More...
 
virtual const
ReferenceTrajectoryCollection 
trajectories (const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const ExternalPredictionCollection &external, const reco::BeamSpot &beamSpot) const
 
virtual ~DualKalmanFactory ()
 
- Public Member Functions inherited from TrajectoryFactoryBase
MaterialEffects materialEffects (void) const
 
PropagationDirection propagationDirection (void) const
 
 TrajectoryFactoryBase (const edm::ParameterSet &config)
 
virtual ~TrajectoryFactoryBase (void)
 

Protected Member Functions

const DualKalmanInput referenceStateAndRecHits (const ConstTrajTrackPair &track) const
 
- Protected Member Functions inherited from TrajectoryFactoryBase
virtual const TrajectoryInput innermostStateAndRecHits (const ConstTrajTrackPair &track) const
 
virtual const
Trajectory::DataContainer 
orderedTrajectoryMeasurements (const Trajectory &trajectory) const
 
bool sameSurface (const Surface &s1, const Surface &s2) const
 
bool useRecHit (const TransientTrackingRecHit::ConstRecHitPointer &hitPtr) const
 

Protected Attributes

const double theMass
 
const int theResidMethod
 
- Protected Attributes inherited from TrajectoryFactoryBase
bool theUseBeamSpot
 

Additional Inherited Members

- Public Types inherited from TrajectoryFactoryBase
typedef
AlignmentAlgorithmBase::ConstTrajTrackPair 
ConstTrajTrackPair
 
typedef
AlignmentAlgorithmBase::ConstTrajTrackPairCollection 
ConstTrajTrackPairCollection
 
typedef std::vector
< TrajectoryStateOnSurface
ExternalPredictionCollection
 
typedef
ReferenceTrajectoryBase::MaterialEffects 
MaterialEffects
 
typedef std::vector
< ReferenceTrajectoryPtr
ReferenceTrajectoryCollection
 
typedef
ReferenceTrajectoryBase::ReferenceTrajectoryPtr 
ReferenceTrajectoryPtr
 
typedef std::pair
< TrajectoryStateOnSurface,
TransientTrackingRecHit::ConstRecHitContainer
TrajectoryInput
 

Detailed Description

A factory that produces reference trajectory instances of class DualKalmanTrajectory from a given TrajTrackPairCollection.

Currently two methods to set residual and error can be choosen via cfg: 1: the unbiased residal approach 2: the pull approach

Author
: Gero Flucke date : October 2008
Revision:
1.4
Date:
2011/12/19 14:40:01

(last update by

Author:
mussgill

)

Definition at line 32 of file DualKalmanFactory.cc.

Constructor & Destructor Documentation

DualKalmanFactory::DualKalmanFactory ( const edm::ParameterSet config)

Definition at line 75 of file DualKalmanFactory.cc.

Referenced by clone().

76  : TrajectoryFactoryBase(config), theMass(config.getParameter<double>("ParticleMass")),
77  theResidMethod(config.getParameter<int>("ResidualMethod"))
78 {
79  // Since theResidMethod is passed to DualKalmanTrajectory, valid values are checked there.
80  edm::LogInfo("Alignment") << "@SUB=DualKalmanFactory" << "Factory created.";
81 }
T getParameter(std::string const &) const
TrajectoryFactoryBase(const edm::ParameterSet &config)
DualKalmanFactory::~DualKalmanFactory ( )
virtual

Definition at line 85 of file DualKalmanFactory.cc.

85 {}

Member Function Documentation

virtual DualKalmanFactory* DualKalmanFactory::clone ( void  ) const
inlinevirtual

Implements TrajectoryFactoryBase.

Definition at line 49 of file DualKalmanFactory.cc.

References DualKalmanFactory().

49 { return new DualKalmanFactory(*this); }
DualKalmanFactory(const edm::ParameterSet &config)
const DualKalmanFactory::DualKalmanInput DualKalmanFactory::referenceStateAndRecHits ( const ConstTrajTrackPair track) const
protected

Definition at line 139 of file DualKalmanFactory.cc.

References DualKalmanFactory::DualKalmanInput::bwdRecHitNums, DualKalmanFactory::DualKalmanInput::fwdRecHitNums, LaserDQM_cfg::input, TrajectoryFactoryBase::orderedTrajectoryMeasurements(), DualKalmanFactory::DualKalmanInput::refTsos, DualKalmanFactory::DualKalmanInput::trajMeasurements, and TrajectoryFactoryBase::useRecHit().

Referenced by trajectories().

140 {
141  // Same idea as in DualTrajectoryFactory::referenceStateAndRecHits(..):
142  // Split trajectory in the middle, take middle as reference and provide first
143  // and second half of hits, each starting from this middle hit.
144  // In contrast to DualTrajectoryFactory we deal here with indices and not rechits directly
145  // to be able to get measurements and uncertainties later from the Trajectory that is
146  // provided by the Kalman track fit.
147 
148  DualKalmanInput input;
149 
150  // get the trajectory measurements in the correct order, i.e. reverse if needed
151  input.trajMeasurements = this->orderedTrajectoryMeasurements(*track.first);
152 
153  // get indices of relevant trajectory measurements to find middle of them
154  std::vector<unsigned int> usedTrajMeasNums;
155  for (unsigned int iM = 0; iM < input.trajMeasurements.size(); ++iM) {
156  if (this->useRecHit(input.trajMeasurements[iM].recHit())) usedTrajMeasNums.push_back(iM);
157  }
158  unsigned int nRefStateMeas = usedTrajMeasNums.size()/2;
159 
160  // get the valid RecHits numbers
161  for (unsigned int iMeas = 0; iMeas < usedTrajMeasNums.size(); ++iMeas) {
162  if (iMeas < nRefStateMeas) {
163  input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
164  } else if (iMeas > nRefStateMeas) {
165  input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
166  } else { // iMeas == nRefStateMeas
167  if (input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState().isValid()) {
168  input.refTsos = input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState();
169  input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
170  input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
171  } else {
172  // if the hit/tsos of the middle hit is not valid, try the next one...
173  ++nRefStateMeas; // but keep hit if only TSOS bad
174  input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
175  }
176  }
177  }
178 
179  // bring input.fwdRecHits into correct order
180  std::reverse(input.bwdRecHitNums.begin(), input.bwdRecHitNums.end());
181 
182  return input;
183 }
virtual const Trajectory::DataContainer orderedTrajectoryMeasurements(const Trajectory &trajectory) const
bool useRecHit(const TransientTrackingRecHit::ConstRecHitPointer &hitPtr) const
const DualKalmanFactory::ReferenceTrajectoryCollection DualKalmanFactory::trajectories ( const edm::EventSetup setup,
const ConstTrajTrackPairCollection tracks,
const reco::BeamSpot beamSpot 
) const
virtual

Produce the reference trajectories.

Implements TrajectoryFactoryBase.

Definition at line 90 of file DualKalmanFactory.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, DualKalmanFactory::DualKalmanInput::bwdRecHitNums, DualKalmanFactory::DualKalmanInput::fwdRecHitNums, edm::EventSetup::get(), LaserDQM_cfg::input, TrajectoryStateOnSurface::isValid(), TrajectoryFactoryBase::materialEffects(), edm::ESHandle< class >::product(), TrajectoryFactoryBase::propagationDirection(), referenceStateAndRecHits(), DualKalmanFactory::DualKalmanInput::refTsos, theMass, theResidMethod, TrajectoryFactoryBase::theUseBeamSpot, and DualKalmanFactory::DualKalmanInput::trajMeasurements.

Referenced by trajectories().

93 {
95 
96  edm::ESHandle<MagneticField> magneticField;
97  setup.get<IdealMagneticFieldRecord>().get(magneticField);
98 
99  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
100 
101  while (itTracks != tracks.end()) {
102  const DualKalmanInput input = this->referenceStateAndRecHits(*itTracks);
103  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
104  if (input.refTsos.isValid()) {
105  ReferenceTrajectoryPtr ptr(new DualKalmanTrajectory(input.trajMeasurements,
106  input.refTsos,
107  input.fwdRecHitNums,
108  input.bwdRecHitNums,
109  magneticField.product(),
110  this->materialEffects(),
111  this->propagationDirection(),
113  theResidMethod));
114  trajectories.push_back(ptr);
115  }
116  ++itTracks;
117  }
118 
119  return trajectories;
120 }
MaterialEffects materialEffects(void) const
virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const
Produce the reference trajectories.
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
PropagationDirection propagationDirection(void) const
const DualKalmanInput referenceStateAndRecHits(const ConstTrajTrackPair &track) const
ReferenceTrajectoryBase::ReferenceTrajectoryPtr ReferenceTrajectoryPtr
const DualKalmanFactory::ReferenceTrajectoryCollection DualKalmanFactory::trajectories ( const edm::EventSetup setup,
const ConstTrajTrackPairCollection tracks,
const ExternalPredictionCollection external,
const reco::BeamSpot beamSpot 
) const
virtual

Implements TrajectoryFactoryBase.

Definition at line 124 of file DualKalmanFactory.cc.

References trajectories().

128 {
130 
131  edm::LogError("Alignment") << "@SUB=DualKalmanFactory::trajectories"
132  << "Not implemented with ExternalPrediction.";
133  return trajectories;
134 }
virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const
Produce the reference trajectories.
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection

Member Data Documentation

const double DualKalmanFactory::theMass
protected

Definition at line 67 of file DualKalmanFactory.cc.

Referenced by twikiExport.Constituent::__str__(), and trajectories().

const int DualKalmanFactory::theResidMethod
protected

Definition at line 68 of file DualKalmanFactory.cc.

Referenced by trajectories().