Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-28 22:48:48

0001 #ifndef RecoParticleFlow_PFTracking_PFDisplacedVertexCandidateFinder_h
0002 #define RecoParticleFlow_PFTracking_PFDisplacedVertexCandidateFinder_h
0003 
0004 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertexCandidate.h"
0005 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertexCandidateFwd.h"
0006 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0007 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0008 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0009 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0010 #include "CommonTools/Utils/interface/KinematicTables.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 
0015 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
0016 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
0017 
0018 /// \brief Displaced Vertex Candidate Finder
0019 /*!
0020   \author Maxime Gouzevitch
0021   \date October 2009
0022 */
0023 
0024 class MagneticField;
0025 
0026 class PFDisplacedVertexCandidateFinder {
0027 public:
0028   PFDisplacedVertexCandidateFinder();
0029 
0030   ~PFDisplacedVertexCandidateFinder();
0031 
0032   /// Mask used to spot if a track is free or not
0033   typedef std::vector<bool> Mask;
0034 
0035   typedef std::list<reco::TrackBaseRef>::iterator IE;
0036   typedef std::list<reco::TrackBaseRef>::const_iterator IEC;
0037   typedef reco::PFDisplacedVertexCandidateCollection::const_iterator IBC;
0038 
0039   /// --------- Set different algo parameters ------ ///
0040 
0041   /// Sets algo parameters for the vertex candidate finder
0042   void setParameters(double dcaCut, double primaryVertexCut, double dcaPInnerHitCut, const edm::ParameterSet& ps_trk) {
0043     dcaCut_ = dcaCut;
0044     primaryVertexCut2_ = primaryVertexCut * primaryVertexCut;
0045     dcaPInnerHitCut2_ = dcaPInnerHitCut * dcaPInnerHitCut;
0046     nChi2_max_ = ps_trk.getParameter<double>("nChi2_max");
0047     pt_min_ = ps_trk.getParameter<double>("pt_min");
0048     pt_min_prim_ = ps_trk.getParameter<double>("pt_min_prim");
0049     dxy_ = ps_trk.getParameter<double>("dxy");
0050     qoverpError_max_ = ps_trk.getParameter<double>("qoverpError_max");
0051   }
0052 
0053   /// sets debug printout flag
0054   void setDebug(bool debug) { debug_ = debug; }
0055 
0056   /// Set the imput collection of tracks and calculate their
0057   /// trajectory parameters the Global Trajectory Parameters
0058   void setInput(const edm::Handle<reco::TrackCollection>& trackh, const MagneticField* magField);
0059 
0060   /// \return auto_ptr to collection of DisplacedVertexCandidates
0061   std::unique_ptr<reco::PFDisplacedVertexCandidateCollection> transferVertexCandidates() {
0062     return std::move(vertexCandidates_);
0063   }
0064 
0065   const std::unique_ptr<reco::PFDisplacedVertexCandidateCollection>& vertexCandidates() const {
0066     return std::move(vertexCandidates_);
0067   }
0068 
0069   /// -------- Main function which find vertices -------- ///
0070 
0071   void findDisplacedVertexCandidates();
0072 
0073   void setPrimaryVertex(edm::Handle<reco::VertexCollection> mainVertexHandle,
0074                         edm::Handle<reco::BeamSpot> beamSpotHandle);
0075 
0076 private:
0077   /// -------- Different steps of the finder algorithm -------- ///
0078 
0079   /// Recursive procedure to associate tracks together
0080   IE associate(IE next, IE last, reco::PFDisplacedVertexCandidate& tempVertexCandidate);
0081 
0082   /// Check whether 2 elements are linked and fill the link parameters
0083   void link(const reco::TrackBaseRef& el1,
0084             const reco::TrackBaseRef& el2,
0085             double& dist,
0086             GlobalPoint& crossing_point,
0087             reco::PFDisplacedVertexCandidate::VertexLinkTest& linktest);
0088 
0089   /// Compute missing links in the displacedVertexCandidates
0090   /// (the recursive procedure does not build all links)
0091   void packLinks(reco::PFDisplacedVertexCandidate& vertexCandidate);
0092 
0093   /// -------- TOOLS -------- //
0094 
0095   /// Allows to calculate the helix aproximation for a given track
0096   /// which may be then extrapolated to any point.
0097   GlobalTrajectoryParameters getGlobalTrajectoryParameters(const reco::Track*) const;
0098 
0099   /// Quality Criterion on the Pt resolution to select a Track
0100   bool goodPtResolution(const reco::TrackBaseRef& trackref) const;
0101 
0102   /// A function which gather the information
0103   /// if a track is available for vertexing
0104   bool isSelected(const reco::TrackBaseRef& trackref) { return goodPtResolution(trackref); }
0105 
0106   friend std::ostream& operator<<(std::ostream&, const PFDisplacedVertexCandidateFinder&);
0107 
0108   /// -------- Members -------- ///
0109 
0110   std::unique_ptr<reco::PFDisplacedVertexCandidateCollection> vertexCandidates_;
0111 
0112   /// The track refs
0113   std::list<reco::TrackBaseRef> eventTracks_;
0114 
0115   /// The trackMask allows to keep the information on the
0116   /// tracks which are still free and those which are already
0117   /// used or disabled.
0118   Mask trackMask_;
0119   /// The Trajectories vector allow to calculate snd to store
0120   /// only once the track trajectory parameters
0121   std::vector<GlobalTrajectoryParameters> eventTrackTrajectories_;
0122 
0123   /// ----- Algo parameters for the vertex finder ---- ///
0124 
0125   /// Distance of minimal approach below which
0126   /// two tracks are considered as linked together
0127   double dcaCut_;
0128   /// Do not reconstruct vertices wich are too close to the beam pipe
0129   double primaryVertexCut2_;
0130   /// Maximum distance between the DCA Point and the inner hit of the track
0131   double dcaPInnerHitCut2_;
0132 
0133   /// Tracks preselection to reduce the combinatorics in PFDisplacedVertexCandidates
0134   /// this cuts are repeated then in a smarter way in the PFDisplacedVertexFinder
0135   /// be sure you are consistent between them
0136   double nChi2_max_;
0137   double pt_min_;
0138 
0139   double pt_min_prim_;
0140   double dxy_;
0141   double qoverpError_max_;
0142   /// Max number of expected vertexCandidates in the event
0143   /// Used to allocate the memory and avoid multiple copy
0144   unsigned vertexCandidatesSize_;
0145 
0146   // Two track minimum distance algo
0147   TwoTrackMinimumDistance theMinimum_;
0148 
0149   math::XYZPoint pvtx_;
0150 
0151   /// if true, debug printouts activated
0152   bool debug_;
0153 
0154   // Tracker geometry for extrapolation
0155   const MagneticField* magField_;
0156 
0157   edm::soa::PtEtaPhiTable track_table_;
0158 };
0159 
0160 #endif