Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:58

0001 /** \class DAFTrackProducerAlgorithm
0002   *  All is needed to run the deterministic annealing algorithm. Ported from ORCA. 
0003   *
0004   *  \author tropiano, genta
0005   *  \review in May 2014 by brondolin 
0006   */
0007 
0008 #ifndef DAFTrackProducerAlgorithm_h
0009 #define DAFTrackProducerAlgorithm_h
0010 
0011 #include "AlgoProductTraits.h"
0012 
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0015 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
0016 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
0017 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0018 #include "RecoTracker/MeasurementDet/interface/MeasurementTrackerEvent.h"
0019 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
0020 
0021 class MagneticField;
0022 class TrackingGeometry;
0023 class TrajAnnealing;
0024 class TrajectoryFitter;
0025 class Trajectory;
0026 class TrajectoryStateOnSurface;
0027 class TransientTrackingRecHitBuilder;
0028 class MultiRecHitCollector;
0029 class SiTrackerMultiRecHitUpdator;
0030 namespace reco {
0031   class Track;
0032 }
0033 
0034 class DAFTrackProducerAlgorithm : public AlgoProductTraits<reco::Track> {
0035 public:
0036   using Base = AlgoProductTraits<reco::Track>;
0037   using TrackCollection = typename Base::TrackCollection;
0038   using AlgoProductCollection = typename Base::AlgoProductCollection;
0039 
0040   using TrajAnnealingCollection = std::vector<TrajAnnealing>;
0041 
0042 public:
0043   DAFTrackProducerAlgorithm(const edm::ParameterSet& conf);
0044   ~DAFTrackProducerAlgorithm() {}
0045 
0046   /// Run the Final Fit taking TrackCandidates as input
0047   void runWithCandidate(const TrackingGeometry*,
0048                         const MagneticField*,
0049                         //const TrackCandidateCollection&,
0050                         const TrajTrackAssociationCollection&,
0051                         const MeasurementTrackerEvent* measTk,
0052                         const TrajectoryFitter*,
0053                         const TransientTrackingRecHitBuilder*,
0054                         const MultiRecHitCollector* measurementTracker,
0055                         const SiTrackerMultiRecHitUpdator*,
0056                         const reco::BeamSpot&,
0057                         AlgoProductCollection&,
0058                         TrajAnnealingCollection&,
0059                         bool,
0060                         AlgoProductCollection&,
0061                         AlgoProductCollection&) const;
0062 
0063 private:
0064   /// Construct Tracks to be put in the event
0065   bool buildTrack(
0066       const Trajectory, AlgoProductCollection& algoResults, float, const reco::BeamSpot&, const reco::TrackRef*) const;
0067 
0068   /// accomplishes the fitting-smoothing step for each annealing value
0069   Trajectory fit(const std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>& hits,
0070                  const TrajectoryFitter* theFitter,
0071                  Trajectory vtraj) const;
0072 
0073   //calculates the ndof according to the DAF prescription
0074   float calculateNdof(const Trajectory vtraj) const;
0075 
0076   //creates MultiRecHits out of a KF trajectory
0077   std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> collectHits(
0078       const Trajectory vtraj,
0079       const MultiRecHitCollector* measurementCollector,
0080       const MeasurementTrackerEvent* measTk) const;
0081 
0082   //updates the hits with the specified annealing factor
0083   std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> updateHits(
0084       const Trajectory vtraj,
0085       const SiTrackerMultiRecHitUpdator* updator,
0086       const MeasurementTrackerEvent* theMTE,
0087       double annealing) const;
0088 
0089   //removes from the trajectory isolated hits with very low weight
0090   void filter(const TrajectoryFitter* fitter,
0091               std::vector<Trajectory>& input,
0092               int minhits,
0093               std::vector<Trajectory>& output,
0094               const TransientTrackingRecHitBuilder* builder) const;
0095 
0096   int countingGoodHits(const Trajectory traj) const;
0097 
0098   int checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const;
0099 
0100   void PrintHit(const TrackingRecHit* const& hit, TrajectoryStateOnSurface& tsos) const;
0101 
0102   edm::ParameterSet conf_;
0103   int minHits_;
0104 };
0105 
0106 #endif