CMS 3D CMS Logo

DAFTrackProducerAlgorithm.h
Go to the documentation of this file.
1 
8 #ifndef DAFTrackProducerAlgorithm_h
9 #define DAFTrackProducerAlgorithm_h
10 
11 #include "AlgoProductTraits.h"
12 
20 
21 class MagneticField;
22 class TrackingGeometry;
23 class TrajAnnealing;
24 class TrajectoryFitter;
25 class Trajectory;
30 namespace reco {
31  class Track;
32 }
33 
34 class DAFTrackProducerAlgorithm : public AlgoProductTraits<reco::Track> {
35 public:
39 
40  using TrajAnnealingCollection = std::vector<TrajAnnealing>;
41 
42 public:
45 
48  const MagneticField*,
49  //const TrackCandidateCollection&,
51  const MeasurementTrackerEvent* measTk,
52  const TrajectoryFitter*,
56  const reco::BeamSpot&,
59  bool,
61  AlgoProductCollection&) const;
62 
63 private:
65  bool buildTrack(
66  const Trajectory, AlgoProductCollection& algoResults, float, const reco::BeamSpot&, const reco::TrackRef*) const;
67 
69  Trajectory fit(const std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>& hits,
70  const TrajectoryFitter* theFitter,
71  Trajectory vtraj) const;
72 
73  //calculates the ndof according to the DAF prescription
74  float calculateNdof(const Trajectory vtraj) const;
75 
76  //creates MultiRecHits out of a KF trajectory
77  std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> collectHits(
78  const Trajectory vtraj,
79  const MultiRecHitCollector* measurementCollector,
80  const MeasurementTrackerEvent* measTk) const;
81 
82  //updates the hits with the specified annealing factor
83  std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> updateHits(
84  const Trajectory vtraj,
86  const MeasurementTrackerEvent* theMTE,
87  double annealing) const;
88 
89  //removes from the trajectory isolated hits with very low weight
90  void filter(const TrajectoryFitter* fitter,
91  std::vector<Trajectory>& input,
92  int minhits,
93  std::vector<Trajectory>& output,
94  const TransientTrackingRecHitBuilder* builder) const;
95 
96  int countingGoodHits(const Trajectory traj) const;
97 
98  int checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const;
99 
100  void PrintHit(const TrackingRecHit* const& hit, TrajectoryStateOnSurface& tsos) const;
101 
103  int minHits_;
104 };
105 
106 #endif
Trajectory fit(const std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > &hits, const TrajectoryFitter *theFitter, Trajectory vtraj) const
accomplishes the fitting-smoothing step for each annealing value
void PrintHit(const TrackingRecHit *const &hit, TrajectoryStateOnSurface &tsos) const
float calculateNdof(const Trajectory vtraj) const
bool buildTrack(const Trajectory, AlgoProductCollection &algoResults, float, const reco::BeamSpot &, const reco::TrackRef *) const
Construct Tracks to be put in the event.
std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > updateHits(const Trajectory vtraj, const SiTrackerMultiRecHitUpdator *updator, const MeasurementTrackerEvent *theMTE, double annealing) const
typename Base::TrackCollection TrackCollection
std::vector< TrajAnnealing > TrajAnnealingCollection
Definition: TrajAnnealing.h:36
std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > collectHits(const Trajectory vtraj, const MultiRecHitCollector *measurementCollector, const MeasurementTrackerEvent *measTk) const
static std::string const input
Definition: EdmProvDump.cc:50
std::vector< TrajAnnealing > TrajAnnealingCollection
std::vector< reco::Track > TrackCollection
void filter(const TrajectoryFitter *fitter, std::vector< Trajectory > &input, int minhits, std::vector< Trajectory > &output, const TransientTrackingRecHitBuilder *builder) const
std::vector< AlgoProduct > AlgoProductCollection
fixed size matrix
int checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const
int countingGoodHits(const Trajectory traj) const
typename Base::AlgoProductCollection AlgoProductCollection
DAFTrackProducerAlgorithm(const edm::ParameterSet &conf)
void runWithCandidate(const TrackingGeometry *, const MagneticField *, const TrajTrackAssociationCollection &, const MeasurementTrackerEvent *measTk, const TrajectoryFitter *, const TransientTrackingRecHitBuilder *, const MultiRecHitCollector *measurementTracker, const SiTrackerMultiRecHitUpdator *, const reco::BeamSpot &, AlgoProductCollection &, TrajAnnealingCollection &, bool, AlgoProductCollection &, AlgoProductCollection &) const
Run the Final Fit taking TrackCandidates as input.