CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DAFTrackProducerAlgorithm.h
Go to the documentation of this file.
1 
8 #ifndef DAFTrackProducerAlgorithm_h
9 #define DAFTrackProducerAlgorithm_h
10 
16 
17 class MagneticField;
18 class TrackingGeometry;
19 class TrajAnnealing;
20 class TrajectoryFitter;
21 class Trajectory;
27 namespace reco{
28  class Track;
29 }
30 
32 
33  typedef std::pair<Trajectory*, std::pair<reco::Track*,PropagationDirection> > AlgoProduct;
34  typedef std::vector< AlgoProduct > AlgoProductCollection;
35  typedef std::vector<TrajAnnealing> TrajAnnealingCollection;
36 
37  public:
38 
41 
43  void runWithCandidate(const TrackingGeometry *,
44  const MagneticField *,
45  //const TrackCandidateCollection&,
46  const std::vector<Trajectory>&,
47  const MeasurementTrackerEvent *measTk,
48  const TrajectoryFitter *,
52  const reco::BeamSpot&,
55  bool ) const;
56 
57  private:
59  bool buildTrack(const Trajectory,
60  AlgoProductCollection& algoResults,
61  float,
62  const reco::BeamSpot&) const;
63 
67  const TrajectoryFitter * theFitter,
68  Trajectory vtraj) const;
69 
70  //calculates the ndof according to the DAF prescription
71  float calculateNdof(const Trajectory vtraj) const;
72 
73  //creates MultiRecHits out of a KF trajectory
74  std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> collectHits(
75  const Trajectory vtraj,
76  const MultiRecHitCollector* measurementCollector,
77  const MeasurementTrackerEvent *measTk ) const;
78 
79  //updates the hits with the specified annealing factor
80  std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> updateHits(
81  const Trajectory vtraj,
83  double annealing) const;
84 
85  //removes from the trajectory isolated hits with very low weight
86  void filter(const TrajectoryFitter* fitter,
87  std::vector<Trajectory>& input,
88  int minhits, std::vector<Trajectory>& output,
89  const TransientTrackingRecHitBuilder* builder) const;
90 
91  int checkHits( Trajectory iInitTraj, const Trajectory iFinalTraj) const;
92 
94  int minHits_;
95 };
96 
97 #endif
std::pair< Trajectory *, std::pair< reco::Track *, PropagationDirection > > AlgoProduct
std::vector< ConstRecHitPointer > RecHitContainer
int checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const
std::vector< TrajAnnealing > TrajAnnealingCollection
Definition: TrajAnnealing.h:42
std::vector< TrajAnnealing > TrajAnnealingCollection
static std::string const input
Definition: EdmProvDump.cc:43
void filter(const TrajectoryFitter *fitter, std::vector< Trajectory > &input, int minhits, std::vector< Trajectory > &output, const TransientTrackingRecHitBuilder *builder) const
std::vector< AlgoProduct > AlgoProductCollection
std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > updateHits(const Trajectory vtraj, const SiTrackerMultiRecHitUpdator *updator, double annealing) const
Trajectory fit(const std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > &hits, const TrajectoryFitter *theFitter, Trajectory vtraj) const
accomplishes the fitting-smoothing step for each annealing value
tuple conf
Definition: dbtoconf.py:185
std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > collectHits(const Trajectory vtraj, const MultiRecHitCollector *measurementCollector, const MeasurementTrackerEvent *measTk) const
bool buildTrack(const Trajectory, AlgoProductCollection &algoResults, float, const reco::BeamSpot &) const
Construct Tracks to be put in the event.
float calculateNdof(const Trajectory vtraj) const
void runWithCandidate(const TrackingGeometry *, const MagneticField *, const std::vector< Trajectory > &, const MeasurementTrackerEvent *measTk, const TrajectoryFitter *, const TransientTrackingRecHitBuilder *, const MultiRecHitCollector *measurementTracker, const SiTrackerMultiRecHitUpdator *, const reco::BeamSpot &, AlgoProductCollection &, TrajAnnealingCollection &, bool) const
Run the Final Fit taking TrackCandidates as input.
DAFTrackProducerAlgorithm(const edm::ParameterSet &conf)