Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:42

0001 #ifndef TRACKINGTOOLS_TRANSIENTRACKBUILDER_H
0002 #define TRACKINGTOOLS_TRANSIENTRACKBUILDER_H
0003 
0004 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0005 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0006 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0007 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0008 #include "DataFormats/Common/interface/ValueMap.h"
0009 
0010 /**
0011    * Helper class to build TransientTrack from the persistent Track.
0012    * This is obtained from the eventSetup, as given in the example in the test
0013    * directory.
0014    */
0015 
0016 class TransientTrackBuilder {
0017 public:
0018   TransientTrackBuilder(const MagneticField* field, const edm::ESHandle<GlobalTrackingGeometry>& trackingGeometry)
0019       : theField(field), theTrackingGeometry(trackingGeometry) {}
0020 
0021   reco::TransientTrack build(const reco::Track* p) const;
0022   reco::TransientTrack build(const reco::Track& p) const;
0023   reco::TransientTrack build(const reco::GsfTrack* p) const;
0024   reco::TransientTrack build(const reco::GsfTrack& p) const;
0025 
0026   reco::TransientTrack build(const reco::TrackRef* p) const;
0027   reco::TransientTrack build(const reco::TrackRef& p) const;
0028   reco::TransientTrack build(const reco::GsfTrackRef* p) const;
0029   reco::TransientTrack build(const reco::GsfTrackRef& p) const;
0030 
0031   reco::TransientTrack build(const reco::CandidatePtr* p) const;
0032   reco::TransientTrack build(const reco::CandidatePtr& p) const;
0033 
0034   std::vector<reco::TransientTrack> build(const edm::Handle<reco::TrackCollection>& trkColl) const;
0035   std::vector<reco::TransientTrack> build(const edm::Handle<reco::GsfTrackCollection>& trkColl) const;
0036   std::vector<reco::TransientTrack> build(const edm::Handle<edm::View<reco::Track> >& trkColl) const;
0037 
0038   std::vector<reco::TransientTrack> build(const edm::Handle<reco::TrackCollection>& trkColl,
0039                                           const edm::ValueMap<float>& trackTimes,
0040                                           const edm::ValueMap<float>& trackTimeResos) const;
0041   std::vector<reco::TransientTrack> build(const edm::Handle<reco::GsfTrackCollection>& trkColl,
0042                                           const edm::ValueMap<float>& trackTimes,
0043                                           const edm::ValueMap<float>& trackTimeResos) const;
0044   std::vector<reco::TransientTrack> build(const edm::Handle<edm::View<reco::Track> >& trkColl,
0045                                           const edm::ValueMap<float>& trackTimes,
0046                                           const edm::ValueMap<float>& trackTimeResos) const;
0047 
0048   std::vector<reco::TransientTrack> build(const edm::Handle<reco::TrackCollection>& trkColl,
0049                                           const reco::BeamSpot& beamSpot) const;
0050   std::vector<reco::TransientTrack> build(const edm::Handle<reco::GsfTrackCollection>& trkColl,
0051                                           const reco::BeamSpot& beamSpot) const;
0052   std::vector<reco::TransientTrack> build(const edm::Handle<edm::View<reco::Track> >& trkColl,
0053                                           const reco::BeamSpot& beamSpot) const;
0054 
0055   std::vector<reco::TransientTrack> build(const edm::Handle<reco::TrackCollection>& trkColl,
0056                                           const reco::BeamSpot& beamSpot,
0057                                           const edm::ValueMap<float>& trackTimes,
0058                                           const edm::ValueMap<float>& trackTimeResos) const;
0059   std::vector<reco::TransientTrack> build(const edm::Handle<reco::GsfTrackCollection>& trkColl,
0060                                           const reco::BeamSpot& beamSpot,
0061                                           const edm::ValueMap<float>& trackTimes,
0062                                           const edm::ValueMap<float>& trackTimeResos) const;
0063   std::vector<reco::TransientTrack> build(const edm::Handle<edm::View<reco::Track> >& trkColl,
0064                                           const reco::BeamSpot& beamSpot,
0065                                           const edm::ValueMap<float>& trackTimes,
0066                                           const edm::ValueMap<float>& trackTimeResos) const;
0067 
0068   reco::TransientTrack build(const FreeTrajectoryState& fts) const;
0069 
0070   const MagneticField* field() const { return theField; }
0071   const edm::ESHandle<GlobalTrackingGeometry> trackingGeometry() const { return theTrackingGeometry; }
0072   static constexpr float defaultInvalidTrackTimeReso = 0.350f;
0073 
0074 private:
0075   const MagneticField* theField;
0076   edm::ESHandle<GlobalTrackingGeometry> theTrackingGeometry;
0077 };
0078 
0079 #endif