Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:28:18

0001 #ifndef DAClusterizerInZ_h
0002 #define DAClusterizerInZ_h
0003 
0004 /**\class DAClusterizerInZ 
0005  
0006   Description: separates event tracks into clusters along the beam line
0007 
0008 */
0009 
0010 #include "RecoVertex/PrimaryVertexProducer/interface/TrackClusterizerInZ.h"
0011 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include <vector>
0014 #include "DataFormats/Math/interface/Error.h"
0015 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0016 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0017 
0018 class DAClusterizerInZ : public TrackClusterizerInZ {
0019 public:
0020   struct track_t {
0021     double z;                        // z-coordinate at point of closest approach to the beamline
0022     double dz2;                      // square of the error of z(pca)
0023     const reco::TransientTrack *tt;  // a pointer to the Transient Track
0024     double Z;                        // Z[i]   for DA clustering
0025     double pi;                       // track weight
0026   };
0027 
0028   struct vertex_t {
0029     double z;   //           z coordinate
0030     double pk;  //           vertex weight for "constrained" clustering
0031     // --- temporary numbers, used during update
0032     double ei;
0033     double sw;
0034     double swz;
0035     double se;
0036     // ---for Tc
0037     double swE;
0038     double Tc;
0039   };
0040 
0041   DAClusterizerInZ(const edm::ParameterSet &conf);
0042 
0043   std::vector<std::vector<reco::TransientTrack> > clusterize(
0044       const std::vector<reco::TransientTrack> &tracks) const override;
0045 
0046   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> &tracks, const int verbosity = 0) const;
0047 
0048   std::vector<track_t> fill(const std::vector<reco::TransientTrack> &tracks) const;
0049 
0050   bool split(double beta, std::vector<track_t> &tks, std::vector<vertex_t> &y, double threshold) const;
0051 
0052   double update(double beta, std::vector<track_t> &tks, std::vector<vertex_t> &y) const;
0053 
0054   double update(double beta, std::vector<track_t> &tks, std::vector<vertex_t> &y, double &) const;
0055 
0056   void dump(const double beta,
0057             const std::vector<vertex_t> &y,
0058             const std::vector<track_t> &tks,
0059             const int verbosity = 0) const;
0060   bool merge(std::vector<vertex_t> &, int) const;
0061   bool merge(std::vector<vertex_t> &, double &) const;
0062   bool purge(std::vector<vertex_t> &, std::vector<track_t> &, double &, const double) const;
0063 
0064   void splitAll(std::vector<vertex_t> &y) const;
0065 
0066   double beta0(const double betamax, std::vector<track_t> &tks, std::vector<vertex_t> &y) const;
0067 
0068   double Eik(const track_t &t, const vertex_t &k) const;
0069 
0070 private:
0071   bool verbose_;
0072   bool useTc_;
0073   float vertexSize_;
0074   int maxIterations_;
0075   double coolingFactor_;
0076   float betamax_;
0077   float betastop_;
0078   double dzCutOff_;
0079   double d0CutOff_;
0080 };
0081 
0082 #endif