CMS 3D CMS Logo

PFTau3ProngSummary.h
Go to the documentation of this file.
1 #ifndef DataFormats_TauReco_PFTau3ProngSummary_h
2 #define DataFormats_TauReco_PFTau3ProngSummary_h
3 
4 /* class PFTau3ProngSummary
5  *
6  * Stores information on the 3 prong summary for a fully reconstructed tau lepton
7  *
8  * author: Ian M. Nugent
9  * The idea of the fully reconstructing the tau using a kinematic fit comes from
10  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
11  * code is a result of the continuation of this work by Ian M. Nugent and Vladimir Cherepanov.
12  */
21 #include "TVector3.h"
22 #include "TLorentzVector.h"
23 #include "TMath.h"
24 
25 namespace reco {
27  public:
29  enum { dimension = 3 };
30  enum { covarianceSize = dimension * ( dimension + 1 ) / 2 };
34 
37  PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,TLorentzVector a1,double vertex_chi2,double vertex_ndf);
38  PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,TLorentzVector a1,double vertex_chi2,double vertex_ndf,TVector3 sv,CovMatrix svcov);
39 
40  virtual ~PFTau3ProngSummary(){}
41 
42  PFTau3ProngSummary* clone() const;
43 
44  virtual bool AddSolution(unsigned int solution, const TLorentzVector& tau, const std::vector<TLorentzVector>& daughter_p4,
45  const std::vector<int>& daughter_charge, const std::vector<int>& daughter_PDGID,
46  bool has3ProngSolution, double solutionChi2, double thetaGJsig);
47 
49  // interface for relevant TIP functions;
50  const VertexRef& primaryVertex() const { return TIP_->primaryVertex(); }
51  CovMatrix primaryVertexCov() const { return TIP_->primaryVertexCov(); }
52  bool hasSecondaryVertex() const { return TIP_->hasSecondaryVertex(); }
53  const VertexRef& secondaryVertex() const { return TIP_->secondaryVertex(); }
54  CovMatrix secondaryVertexCov() const { return TIP_->secondaryVertexCov(); }
55  const Vector& flightLength() const { return TIP_->flightLength(); }
56  double flightLengthSig() const { return TIP_->flightLengthSig(); }
57  CovMatrix flightLenghtCov() const { return TIP_->flightLengthCov(); }
58 
59  // Tau 3 prong functions
60  const TLorentzVector& A1_LV()const { return a1_; }
61  double M_A1() const { return a1_.M(); }
62  double M_12() const; //pi-pi-
63  double M_13() const; //pi-pi+ Dalitz masses
64  double M_23() const; //pi-pi+ Dalitz masses
65  int Tau_Charge() const;
66  const TVector3& HelixFitSecondaryVertex() const{return sv_; }
67  const CovMatrix& HelixFitSecondaryVertexCov() const{return svcov_; }
68  double Vertex_chi2() const { return vertex_chi2_; }
69  double Vertex_ndf() const { return vertex_ndf_; }
70  double Vertex_Prob() const { return TMath::Prob(vertex_chi2_, vertex_ndf_); }
71  bool has3ProngSolution(unsigned int i) const { return has3ProngSolution_[i]; }
72  double Solution_Chi2(unsigned int i) const { return solution_Chi2_[i]; }
73  double SignificanceOfThetaGJ(unsigned int i) const { return thetaGJsig_[i]; } // 0 or less means the theta_GF has
74  // a physical solution
75  const TLorentzVector& Tau(unsigned int i) const { return tau_p4_[i]; }
76  const std::vector<int>& Daughter_PDGID(unsigned int i) const { return daughter_PDGID_[i]; }
77  const std::vector<int>& Daughter_Charge(unsigned int i) const { return daughter_charge_[i]; }
78  const std::vector<TLorentzVector>& Daughter_P4(unsigned int i) const { return daughter_p4_[i]; }
79 
80  private:
82  TLorentzVector a1_;
83  TVector3 sv_;
84  CovMatrix svcov_;
85  double vertex_chi2_;
86  double vertex_ndf_;
87  std::vector<bool> has3ProngSolution_;
88  std::vector<double> solution_Chi2_;
89  std::vector<double> thetaGJsig_;
90  std::vector<TLorentzVector> tau_p4_;
91  std::vector<std::vector<int> > daughter_PDGID_;
92  std::vector<std::vector<int> > daughter_charge_;
93  std::vector<std::vector<TLorentzVector> > daughter_p4_;
94 
95  };
96 }
97 
98 #endif
double SignificanceOfThetaGJ(unsigned int i) const
reco::PFTauTransverseImpactParameterRef TIP_
math::Error< dimension >::type CovMatrix
std::vector< bool > has3ProngSolution_
std::vector< double > thetaGJsig_
std::vector< double > solution_Chi2_
ErrorD< N >::type type
Definition: Error.h:33
const Vector & flightLength() const
const reco::PFTauTransverseImpactParameterRef & PFTauTIP() const
const std::vector< int > & Daughter_PDGID(unsigned int i) const
const TLorentzVector & Tau(unsigned int i) const
const TLorentzVector & A1_LV() const
const TVector3 & HelixFitSecondaryVertex() const
bool has3ProngSolution(unsigned int i) const
std::vector< TLorentzVector > tau_p4_
std::vector< std::vector< int > > daughter_charge_
const std::vector< TLorentzVector > & Daughter_P4(unsigned int i) const
PFTau3ProngSummary()
constructor from values
CovMatrix primaryVertexCov() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
virtual bool AddSolution(unsigned int solution, const TLorentzVector &tau, const std::vector< TLorentzVector > &daughter_p4, const std::vector< int > &daughter_charge, const std::vector< int > &daughter_PDGID, bool has3ProngSolution, double solutionChi2, double thetaGJsig)
CovMatrix flightLenghtCov() const
CovMatrix secondaryVertexCov() const
PFTau3ProngSummary * clone() const
const CovMatrix & HelixFitSecondaryVertexCov() const
fixed size matrix
std::vector< std::vector< int > > daughter_PDGID_
std::vector< std::vector< TLorentzVector > > daughter_p4_
double Solution_Chi2(unsigned int i) const
const VertexRef & primaryVertex() const
const std::vector< int > & Daughter_Charge(unsigned int i) const
const VertexRef & secondaryVertex() const
double flightLengthSig() const