Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:57

0001 #ifndef DataFormats_TauReco_PFTau3ProngSummary_h
0002 #define DataFormats_TauReco_PFTau3ProngSummary_h
0003 
0004 /* class PFTau3ProngSummary
0005  * 
0006  * Stores information on the 3 prong summary for a fully reconstructed tau lepton
0007  *
0008  * author: Ian M. Nugent
0009  * The idea of the fully reconstructing the tau using a kinematic fit comes from
0010  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
0011  * code is a result of the continuation of this work by Ian M. Nugent and Vladimir Cherepanov.
0012  */
0013 #include "DataFormats/Math/interface/Vector3D.h"
0014 #include "DataFormats/Math/interface/Point3D.h"
0015 #include "DataFormats/Math/interface/Error.h"
0016 #include "DataFormats/VertexReco/interface/Vertex.h"
0017 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0018 #include "DataFormats/TauReco/interface/PFTauTransverseImpactParameter.h"
0019 #include "DataFormats/TauReco/interface/PFTauTransverseImpactParameterFwd.h"
0020 #include "TVector3.h"
0021 #include "TLorentzVector.h"
0022 #include "TMath.h"
0023 
0024 namespace reco {
0025   class PFTau3ProngSummary {
0026   public:
0027     enum { ambiguity, minus, plus, nsolutions };
0028     enum { dimension = 3 };
0029     enum { covarianceSize = dimension * (dimension + 1) / 2 };
0030     typedef math::Error<dimension>::type CovMatrix;
0031     typedef math::XYZPoint Point;
0032     typedef math::XYZVector Vector;
0033 
0034     /// constructor from values
0035     PFTau3ProngSummary();
0036     PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,
0037                        TLorentzVector a1,
0038                        double vertex_chi2,
0039                        double vertex_ndf);
0040     PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,
0041                        TLorentzVector a1,
0042                        double vertex_chi2,
0043                        double vertex_ndf,
0044                        TVector3 sv,
0045                        CovMatrix svcov);
0046 
0047     virtual ~PFTau3ProngSummary() {}
0048 
0049     PFTau3ProngSummary* clone() const;
0050 
0051     virtual bool AddSolution(unsigned int solution,
0052                              const TLorentzVector& tau,
0053                              const std::vector<TLorentzVector>& daughter_p4,
0054                              const std::vector<int>& daughter_charge,
0055                              const std::vector<int>& daughter_PDGID,
0056                              bool has3ProngSolution,
0057                              double solutionChi2,
0058                              double thetaGJsig);
0059 
0060     const reco::PFTauTransverseImpactParameterRef& PFTauTIP() const { return TIP_; }
0061     // interface for relevant TIP functions;
0062     const VertexRef& primaryVertex() const { return TIP_->primaryVertex(); }
0063     CovMatrix primaryVertexCov() const { return TIP_->primaryVertexCov(); }
0064     bool hasSecondaryVertex() const { return TIP_->hasSecondaryVertex(); }
0065     const VertexRef& secondaryVertex() const { return TIP_->secondaryVertex(); }
0066     CovMatrix secondaryVertexCov() const { return TIP_->secondaryVertexCov(); }
0067     const Vector& flightLength() const { return TIP_->flightLength(); }
0068     double flightLengthSig() const { return TIP_->flightLengthSig(); }
0069     CovMatrix flightLenghtCov() const { return TIP_->flightLengthCov(); }
0070 
0071     // Tau 3 prong functions
0072     const TLorentzVector& A1_LV() const { return a1_; }
0073     double M_A1() const { return a1_.M(); }
0074     double M_12() const;  //pi-pi-
0075     double M_13() const;  //pi-pi+ Dalitz masses
0076     double M_23() const;  //pi-pi+ Dalitz masses
0077     int Tau_Charge() const;
0078     const TVector3& HelixFitSecondaryVertex() const { return sv_; }
0079     const CovMatrix& HelixFitSecondaryVertexCov() const { return svcov_; }
0080     double Vertex_chi2() const { return vertex_chi2_; }
0081     double Vertex_ndf() const { return vertex_ndf_; }
0082     double Vertex_Prob() const { return TMath::Prob(vertex_chi2_, vertex_ndf_); }
0083     bool has3ProngSolution(unsigned int i) const { return has3ProngSolution_[i]; }
0084     double Solution_Chi2(unsigned int i) const { return solution_Chi2_[i]; }
0085     double SignificanceOfThetaGJ(unsigned int i) const { return thetaGJsig_[i]; }  // 0 or less means the theta_GF has
0086     // a physical solution
0087     const TLorentzVector& Tau(unsigned int i) const { return tau_p4_[i]; }
0088     const std::vector<int>& Daughter_PDGID(unsigned int i) const { return daughter_PDGID_[i]; }
0089     const std::vector<int>& Daughter_Charge(unsigned int i) const { return daughter_charge_[i]; }
0090     const std::vector<TLorentzVector>& Daughter_P4(unsigned int i) const { return daughter_p4_[i]; }
0091 
0092   private:
0093     reco::PFTauTransverseImpactParameterRef TIP_;
0094     TLorentzVector a1_;
0095     TVector3 sv_;
0096     CovMatrix svcov_;
0097     double vertex_chi2_;
0098     double vertex_ndf_;
0099     std::vector<bool> has3ProngSolution_;
0100     std::vector<double> solution_Chi2_;
0101     std::vector<double> thetaGJsig_;
0102     std::vector<TLorentzVector> tau_p4_;
0103     std::vector<std::vector<int> > daughter_PDGID_;
0104     std::vector<std::vector<int> > daughter_charge_;
0105     std::vector<std::vector<TLorentzVector> > daughter_p4_;
0106   };
0107 }  // namespace reco
0108 
0109 #endif