Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:23

0001 #ifndef VertexReco_Vertex_h
0002 #define VertexReco_Vertex_h
0003 /** \class reco::Vertex
0004  *  
0005  * A reconstructed Vertex providing position, error, chi2, ndof 
0006  * and reconstrudted tracks.
0007  * The vertex can be valid, fake, or invalid.
0008  * A valid vertex is one which has been obtained from a vertex fit of tracks, 
0009  * and all data is meaningful
0010  * A fake vertex is a vertex which was not made out of a proper fit with
0011  * tracks, but still has a position and error (chi2 and ndof are null). 
0012  * For a primary vertex, it could simply be the beam line.
0013  * A fake vertex is considered valid.
0014  * An invalid vertex has no meaningful data.
0015  *
0016  * \author Luca Lista, INFN
0017  *
0018  *
0019  */
0020 #include <Rtypes.h>
0021 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0022 #include "DataFormats/Math/interface/Error.h"
0023 #include "DataFormats/Math/interface/Point3D.h"
0024 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0025 #include "DataFormats/TrackReco/interface/Track.h"
0026 #include "DataFormats/Common/interface/RefToBase.h"
0027 #include <Math/GenVector/PxPyPzE4D.h>
0028 #include <Math/GenVector/PxPyPzM4D.h>
0029 #include "DataFormats/Math/interface/LorentzVector.h"
0030 
0031 namespace reco {
0032 
0033   class Track;
0034 
0035   class Vertex {
0036   public:
0037     /// The iteratator for the vector<TrackRef>
0038     typedef std::vector<TrackBaseRef>::const_iterator trackRef_iterator;
0039     /// point in the space
0040     typedef math::XYZPoint Point;
0041     /// error matrix dimension
0042     constexpr static int dimension = 3;
0043     constexpr static int dimension4D = 4;
0044     /// covariance error matrix (3x3)
0045     typedef math::Error<dimension>::type Error;
0046     /// covariance error matrix (3x3)
0047     typedef math::Error<dimension>::type CovarianceMatrix;
0048     /// covariance error matrix (4x4)
0049     typedef math::Error<dimension4D>::type Error4D;
0050     /// covariance error matrix (4x4)
0051     typedef math::Error<dimension4D>::type CovarianceMatrix4D;
0052     /// matix size
0053     constexpr static int size = dimension * (dimension + 1) / 2, size4D = (dimension4D) * (dimension4D + 1) / 2;
0054     /// index type
0055     typedef unsigned int index;
0056     /// default constructor - The vertex will not be valid. Position, error,
0057     /// chi2, ndof will have random entries, and the vectors of tracks will be empty
0058     /// Use the isValid method to check that your vertex is valid.
0059     Vertex() : chi2_(0.0), ndof_(0), position_(0., 0., 0.), time_(0.) {
0060       validity_ = false;
0061       for (int i = 0; i < size4D; ++i)
0062         covariance_[i] = 0.;
0063     }
0064     /// Constructor for a fake vertex.
0065     Vertex(const Point &, const Error &);
0066     /// Constructor for a fake vertex. 4D
0067     Vertex(const Point &, const Error4D &, double);
0068     /// constructor for a valid vertex, with all data
0069     Vertex(const Point &, const Error &, double chi2, double ndof, size_t size);
0070     /// constructor for a valid vertex, with all data 4D
0071     Vertex(const Point &, const Error4D &, double time, double chi2, double ndof, size_t size);
0072     /// Tells whether the vertex is valid.
0073     bool isValid() const { return validity_; }
0074     /// Tells whether a Vertex is fake, i.e. not a vertex made out of a proper
0075     /// fit with tracks.
0076     /// For a primary vertex, it could simply be the beam line.
0077     bool isFake() const { return (chi2_ == 0 && ndof_ == 0 && tracks_.empty()); }
0078     /// reserve space for the tracks
0079     void reserve(int size, bool refitAsWell = false) {
0080       tracks_.reserve(size);
0081       if (refitAsWell)
0082         refittedTracks_.reserve(size);
0083       weights_.reserve(size);
0084     }
0085     /// add a reference to a Track
0086     template <typename Ref>
0087     void add(Ref const &r, float w = 1.0) {
0088       tracks_.emplace_back(r);
0089       weights_.emplace_back(w * 255.f);
0090     }
0091     /// add the original a Track(reference) and the smoothed Track
0092     void add(const TrackBaseRef &r, const Track &refTrack, float w = 1.0);
0093     void removeTracks();
0094 
0095     ///returns the weight with which a Track has contributed to the vertex-fit.
0096     template <typename TREF>
0097     float trackWeight(const TREF &r) const {
0098       int i = 0;
0099       for (auto const &t : tracks_) {
0100         if ((r.id() == t.id()) && (t.key() == r.key()))
0101           return weights_[i] / 255.f;
0102         ++i;
0103       }
0104       return 0;
0105     }
0106     // track collections
0107     auto const &tracks() const { return tracks_; }
0108     /// first iterator over tracks
0109     trackRef_iterator tracks_begin() const { return tracks_.begin(); }
0110     /// last iterator over tracks
0111     trackRef_iterator tracks_end() const { return tracks_.end(); }
0112     /// number of tracks
0113     size_t tracksSize() const { return tracks_.size(); }
0114     /// python friendly track getting
0115     const TrackBaseRef &trackRefAt(size_t idx) const { return tracks_[idx]; }
0116     /// chi-squares
0117     double chi2() const { return chi2_; }
0118     /** Number of degrees of freedom
0119      *  Meant to be Double32_t for soft-assignment fitters: 
0120      *  tracks may contribute to the vertex with fractional weights.
0121      *  The ndof is then = to the sum of the track weights.
0122      *  see e.g. CMS NOTE-2006/032, CMS NOTE-2004/002
0123      */
0124     double ndof() const { return ndof_; }
0125     /// chi-squared divided by n.d.o.f.
0126     double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
0127     /// position
0128     const Point &position() const { return position_; }
0129     /// x coordinate
0130     double x() const { return position_.X(); }
0131     /// y coordinate
0132     double y() const { return position_.Y(); }
0133     /// z coordinate
0134     double z() const { return position_.Z(); }
0135     /// t coordinate
0136     double t() const { return time_; }
0137     /// error on x
0138     double xError() const { return sqrt(covariance(0, 0)); }
0139     /// error on y
0140     double yError() const { return sqrt(covariance(1, 1)); }
0141     /// error on z
0142     double zError() const { return sqrt(covariance(2, 2)); }
0143     /// error on t
0144     double tError() const { return sqrt(covariance(3, 3)); }
0145     /// (i, j)-th element of error matrix, i, j = 0, ... 2
0146     // Note that:
0147     //   double error( int i, int j ) const
0148     // is OBSOLETE, use covariance(i, j)
0149     double covariance(int i, int j) const { return covariance_[idx(i, j)]; }
0150     /// return SMatrix
0151     CovarianceMatrix covariance() const {
0152       Error m;
0153       fill(m);
0154       return m;
0155     }
0156     /// return SMatrix 4D
0157     CovarianceMatrix4D covariance4D() const {
0158       Error4D m;
0159       fill(m);
0160       return m;
0161     }
0162 
0163     /// return SMatrix
0164     Error error() const {
0165       Error m;
0166       fill(m);
0167       return m;
0168     }
0169     /// return SMatrix
0170     Error4D error4D() const {
0171       Error4D m;
0172       fill(m);
0173       return m;
0174     }
0175 
0176     /// fill SMatrix
0177     void fill(CovarianceMatrix &v) const;
0178     /// 4D version
0179     void fill(CovarianceMatrix4D &v) const;
0180 
0181     /// Checks whether refitted tracks are stored.
0182     bool hasRefittedTracks() const { return !refittedTracks_.empty(); }
0183 
0184     /// Returns the original track which corresponds to a particular refitted Track
0185     /// Throws an exception if now refitted tracks are stored ot the track is not found in the list
0186     TrackBaseRef originalTrack(const Track &refTrack) const;
0187 
0188     /// Returns the refitted track which corresponds to a particular original Track
0189     /// Throws an exception if now refitted tracks are stored ot the track is not found in the list
0190     Track refittedTrack(const TrackBaseRef &track) const;
0191 
0192     /// Returns the refitted track which corresponds to a particular original Track
0193     /// Throws an exception if now refitted tracks are stored ot the track is not found in the list
0194     Track refittedTrack(const TrackRef &track) const;
0195 
0196     /// Returns the container of refitted tracks
0197     const std::vector<Track> &refittedTracks() const { return refittedTracks_; }
0198 
0199     /// Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products
0200     math::XYZTLorentzVectorD p4(float mass = 0.13957018, float minWeight = 0.5) const;
0201 
0202     /// Returns the number of tracks in the vertex with weight above minWeight
0203     unsigned int nTracks(float minWeight = 0.5) const;
0204 
0205     class TrackEqual {
0206     public:
0207       TrackEqual(const Track &t) : track_(t) {}
0208       bool operator()(const Track &t) const { return t.pt() == track_.pt(); }
0209 
0210     private:
0211       const Track &track_;
0212     };
0213 
0214   private:
0215     /// chi-sqared
0216     float chi2_;
0217     /// number of degrees of freedom
0218     float ndof_;
0219     /// position
0220     Point position_;
0221     /// covariance matrix (4x4) as vector
0222     float covariance_[size4D];
0223     /// reference to tracks
0224     std::vector<TrackBaseRef> tracks_;
0225     /// The vector of refitted tracks
0226     std::vector<Track> refittedTracks_;
0227     std::vector<uint8_t> weights_;
0228     /// tells wether the vertex is really valid.
0229     bool validity_;
0230     double time_;
0231 
0232     /// position index
0233     index idx(index i, index j) const {
0234       int a = (i <= j ? i : j), b = (i <= j ? j : i);
0235       return b * (b + 1) / 2 + a;
0236     }
0237   };
0238 
0239 }  // namespace reco
0240 
0241 #endif