CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  */
20 #include "TVector3.h"
21 #include "TLorentzVector.h"
22 #include "TMath.h"
23 
24 namespace reco {
26  public:
28  enum { dimension = 3 };
29  enum { covarianceSize = dimension * (dimension + 1) / 2 };
33 
37  TLorentzVector a1,
38  double vertex_chi2,
39  double vertex_ndf);
41  TLorentzVector a1,
42  double vertex_chi2,
43  double vertex_ndf,
44  TVector3 sv,
45  CovMatrix svcov);
46 
47  virtual ~PFTau3ProngSummary() {}
48 
49  PFTau3ProngSummary* clone() const;
50 
51  virtual bool AddSolution(unsigned int solution,
52  const TLorentzVector& tau,
53  const std::vector<TLorentzVector>& daughter_p4,
54  const std::vector<int>& daughter_charge,
55  const std::vector<int>& daughter_PDGID,
56  bool has3ProngSolution,
57  double solutionChi2,
58  double thetaGJsig);
59 
61  // interface for relevant TIP functions;
62  const VertexRef& primaryVertex() const { return TIP_->primaryVertex(); }
63  CovMatrix primaryVertexCov() const { return TIP_->primaryVertexCov(); }
64  bool hasSecondaryVertex() const { return TIP_->hasSecondaryVertex(); }
65  const VertexRef& secondaryVertex() const { return TIP_->secondaryVertex(); }
66  CovMatrix secondaryVertexCov() const { return TIP_->secondaryVertexCov(); }
67  const Vector& flightLength() const { return TIP_->flightLength(); }
68  double flightLengthSig() const { return TIP_->flightLengthSig(); }
69  CovMatrix flightLenghtCov() const { return TIP_->flightLengthCov(); }
70 
71  // Tau 3 prong functions
72  const TLorentzVector& A1_LV() const { return a1_; }
73  double M_A1() const { return a1_.M(); }
74  double M_12() const; //pi-pi-
75  double M_13() const; //pi-pi+ Dalitz masses
76  double M_23() const; //pi-pi+ Dalitz masses
77  int Tau_Charge() const;
78  const TVector3& HelixFitSecondaryVertex() const { return sv_; }
79  const CovMatrix& HelixFitSecondaryVertexCov() const { return svcov_; }
80  double Vertex_chi2() const { return vertex_chi2_; }
81  double Vertex_ndf() const { return vertex_ndf_; }
82  double Vertex_Prob() const { return TMath::Prob(vertex_chi2_, vertex_ndf_); }
83  bool has3ProngSolution(unsigned int i) const { return has3ProngSolution_[i]; }
84  double Solution_Chi2(unsigned int i) const { return solution_Chi2_[i]; }
85  double SignificanceOfThetaGJ(unsigned int i) const {
86  return thetaGJsig_[i];
87  } // 0 or less means the theta_GF has
88  // a physical solution
89  const TLorentzVector& Tau(unsigned int i) const { return tau_p4_[i]; }
90  const std::vector<int>& Daughter_PDGID(unsigned int i) const { return daughter_PDGID_[i]; }
91  const std::vector<int>& Daughter_Charge(unsigned int i) const { return daughter_charge_[i]; }
92  const std::vector<TLorentzVector>& Daughter_P4(unsigned int i) const { return daughter_p4_[i]; }
93 
94  private:
96  TLorentzVector a1_;
97  TVector3 sv_;
98  CovMatrix svcov_;
99  double vertex_chi2_;
100  double vertex_ndf_;
101  std::vector<bool> has3ProngSolution_;
102  std::vector<double> solution_Chi2_;
103  std::vector<double> thetaGJsig_;
104  std::vector<TLorentzVector> tau_p4_;
105  std::vector<std::vector<int> > daughter_PDGID_;
106  std::vector<std::vector<int> > daughter_charge_;
107  std::vector<std::vector<TLorentzVector> > daughter_p4_;
108  };
109 } // namespace reco
110 
111 #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:32
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:31
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
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