CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
reco::PFTau3ProngSummary Class Reference

#include <PFTau3ProngSummary.h>

Public Types

enum  { ambiguity, minus, plus, nsolutions }
 
enum  { dimension = 3 }
 
enum  { covarianceSize = dimension * (dimension + 1) / 2 }
 
typedef math::Error< dimension >::type CovMatrix
 
typedef math::XYZPoint Point
 
typedef math::XYZVector Vector
 

Public Member Functions

const TLorentzVector & A1_LV () const
 
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)
 
PFTau3ProngSummaryclone () const
 
const std::vector< int > & Daughter_Charge (unsigned int i) const
 
const std::vector< TLorentzVector > & Daughter_P4 (unsigned int i) const
 
const std::vector< int > & Daughter_PDGID (unsigned int i) const
 
CovMatrix flightLenghtCov () const
 
const VectorflightLength () const
 
double flightLengthSig () const
 
bool has3ProngSolution (unsigned int i) const
 
bool hasSecondaryVertex () const
 
const TVector3 & HelixFitSecondaryVertex () const
 
const CovMatrixHelixFitSecondaryVertexCov () const
 
double M_12 () const
 
double M_13 () const
 
double M_23 () const
 
double M_A1 () const
 
 PFTau3ProngSummary ()
 constructor from values More...
 
 PFTau3ProngSummary (reco::PFTauTransverseImpactParameterRef TIP, TLorentzVector a1, double vertex_chi2, double vertex_ndf)
 
 PFTau3ProngSummary (reco::PFTauTransverseImpactParameterRef TIP, TLorentzVector a1, double vertex_chi2, double vertex_ndf, TVector3 sv, CovMatrix svcov)
 
const reco::PFTauTransverseImpactParameterRefPFTauTIP () const
 
const VertexRefprimaryVertex () const
 
CovMatrix primaryVertexCov () const
 
const VertexRefsecondaryVertex () const
 
CovMatrix secondaryVertexCov () const
 
double SignificanceOfThetaGJ (unsigned int i) const
 
double Solution_Chi2 (unsigned int i) const
 
const TLorentzVector & Tau (unsigned int i) const
 
int Tau_Charge () const
 
double Vertex_chi2 () const
 
double Vertex_ndf () const
 
double Vertex_Prob () const
 
virtual ~PFTau3ProngSummary ()
 

Private Attributes

TLorentzVector a1_
 
std::vector< std::vector< int > > daughter_charge_
 
std::vector< std::vector< TLorentzVector > > daughter_p4_
 
std::vector< std::vector< int > > daughter_PDGID_
 
std::vector< bool > has3ProngSolution_
 
std::vector< double > solution_Chi2_
 
TVector3 sv_
 
CovMatrix svcov_
 
std::vector< TLorentzVector > tau_p4_
 
std::vector< double > thetaGJsig_
 
reco::PFTauTransverseImpactParameterRef TIP_
 
double vertex_chi2_
 
double vertex_ndf_
 

Detailed Description

Definition at line 26 of file PFTau3ProngSummary.h.

Member Typedef Documentation

Definition at line 31 of file PFTau3ProngSummary.h.

Definition at line 32 of file PFTau3ProngSummary.h.

Definition at line 33 of file PFTau3ProngSummary.h.

Member Enumeration Documentation

anonymous enum
anonymous enum
Enumerator
dimension 

Definition at line 29 of file PFTau3ProngSummary.h.

anonymous enum

Constructor & Destructor Documentation

PFTau3ProngSummary::PFTau3ProngSummary ( )

constructor from values

Definition at line 51 of file PFTau3ProngSummary.cc.

References daughter_charge_, daughter_p4_, daughter_PDGID_, has3ProngSolution_, mps_fire::i, nsolutions, solution_Chi2_, tau_p4_, and thetaGJsig_.

Referenced by clone().

51  {
52  for (unsigned int i = 0; i < nsolutions; i++) {
53  has3ProngSolution_.push_back(false);
54  solution_Chi2_.push_back(0);
55  thetaGJsig_.push_back(0);
56  tau_p4_.push_back(TLorentzVector(0, 0, 0, 0));
57  daughter_PDGID_.push_back(std::vector<int>());
58  daughter_charge_.push_back(std::vector<int>());
59  daughter_p4_.push_back(std::vector<TLorentzVector>());
60  }
61 }
std::vector< bool > has3ProngSolution_
std::vector< double > thetaGJsig_
std::vector< double > solution_Chi2_
std::vector< TLorentzVector > tau_p4_
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< int > > daughter_PDGID_
std::vector< std::vector< TLorentzVector > > daughter_p4_
PFTau3ProngSummary::PFTau3ProngSummary ( reco::PFTauTransverseImpactParameterRef  TIP,
TLorentzVector  a1,
double  vertex_chi2,
double  vertex_ndf 
)

Definition at line 7 of file PFTau3ProngSummary.cc.

References a1_, daughter_charge_, daughter_p4_, daughter_PDGID_, has3ProngSolution_, mps_fire::i, nsolutions, solution_Chi2_, sv_, svcov_, tau_p4_, thetaGJsig_, TIP_, vertex_chi2_, and vertex_ndf_.

10  {
11  TIP_ = TIP;
12  for (unsigned int i = 0; i < nsolutions; i++) {
13  has3ProngSolution_.push_back(false);
14  solution_Chi2_.push_back(0);
15  thetaGJsig_.push_back(0);
16  tau_p4_.push_back(TLorentzVector(0, 0, 0, 0));
17  daughter_PDGID_.push_back(std::vector<int>());
18  daughter_charge_.push_back(std::vector<int>());
19  daughter_p4_.push_back(std::vector<TLorentzVector>());
20  }
21  a1_ = a1;
22  sv_ = TVector3(TIP_->secondaryVertex()->x(), TIP_->secondaryVertex()->y(), TIP_->secondaryVertex()->z());
23  svcov_ = TIP_->secondaryVertexCov();
24  vertex_chi2_ = vertex_chi2;
25  vertex_ndf_ = vertex_ndf;
26 }
reco::PFTauTransverseImpactParameterRef TIP_
std::vector< bool > has3ProngSolution_
std::vector< double > thetaGJsig_
std::vector< double > solution_Chi2_
std::vector< TLorentzVector > tau_p4_
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< int > > daughter_PDGID_
std::vector< std::vector< TLorentzVector > > daughter_p4_
PFTau3ProngSummary::PFTau3ProngSummary ( reco::PFTauTransverseImpactParameterRef  TIP,
TLorentzVector  a1,
double  vertex_chi2,
double  vertex_ndf,
TVector3  sv,
CovMatrix  svcov 
)

Definition at line 28 of file PFTau3ProngSummary.cc.

References a1_, daughter_charge_, daughter_p4_, daughter_PDGID_, has3ProngSolution_, mps_fire::i, nsolutions, solution_Chi2_, pfDeepBoostedJetPreprocessParams_cfi::sv, sv_, svcov_, tau_p4_, thetaGJsig_, TIP_, vertex_chi2_, and vertex_ndf_.

33  {
34  TIP_ = TIP;
35  for (unsigned int i = 0; i < nsolutions; i++) {
36  has3ProngSolution_.push_back(false);
37  solution_Chi2_.push_back(0);
38  thetaGJsig_.push_back(0);
39  tau_p4_.push_back(TLorentzVector(0, 0, 0, 0));
40  daughter_PDGID_.push_back(std::vector<int>());
41  daughter_charge_.push_back(std::vector<int>());
42  daughter_p4_.push_back(std::vector<TLorentzVector>());
43  }
44  a1_ = a1;
45  sv_ = sv;
46  svcov_ = svcov;
47  vertex_chi2_ = vertex_chi2;
48  vertex_ndf_ = vertex_ndf;
49 }
reco::PFTauTransverseImpactParameterRef TIP_
std::vector< bool > has3ProngSolution_
std::vector< double > thetaGJsig_
std::vector< double > solution_Chi2_
std::vector< TLorentzVector > tau_p4_
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< int > > daughter_PDGID_
std::vector< std::vector< TLorentzVector > > daughter_p4_
virtual reco::PFTau3ProngSummary::~PFTau3ProngSummary ( )
inlinevirtual

Definition at line 48 of file PFTau3ProngSummary.h.

References AddSolution(), clone(), has3ProngSolution(), and metsig::tau.

48 {}

Member Function Documentation

const TLorentzVector& reco::PFTau3ProngSummary::A1_LV ( ) const
inline

Definition at line 73 of file PFTau3ProngSummary.h.

References a1_.

73 { return a1_; }
bool PFTau3ProngSummary::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 
)
virtual

Definition at line 65 of file PFTau3ProngSummary.cc.

References daughter_charge_, daughter_p4_, daughter_PDGID_, has3ProngSolution_, nsolutions, solution_Chi2_, metsig::tau, tau_p4_, and thetaGJsig_.

Referenced by PFTau3ProngReco::produce(), and ~PFTau3ProngSummary().

72  {
73  if (solution < nsolutions) {
74  has3ProngSolution_[solution] = true;
75  solution_Chi2_[solution] = solutionChi2;
76  thetaGJsig_[solution] = thetaGJsig;
77  tau_p4_[solution] = tau;
78  daughter_PDGID_[solution] = daughter_PDGID;
79  daughter_charge_[solution] = daughter_charge;
80  daughter_p4_[solution] = daughter_p4;
81  return true;
82  }
83  return false;
84 }
std::vector< bool > has3ProngSolution_
std::vector< double > thetaGJsig_
std::vector< double > solution_Chi2_
std::vector< TLorentzVector > tau_p4_
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< int > > daughter_PDGID_
std::vector< std::vector< TLorentzVector > > daughter_p4_
PFTau3ProngSummary * PFTau3ProngSummary::clone ( void  ) const

Definition at line 63 of file PFTau3ProngSummary.cc.

References PFTau3ProngSummary().

Referenced by ~PFTau3ProngSummary().

63 { return new PFTau3ProngSummary(*this); }
PFTau3ProngSummary()
constructor from values
const std::vector<int>& reco::PFTau3ProngSummary::Daughter_Charge ( unsigned int  i) const
inline

Definition at line 92 of file PFTau3ProngSummary.h.

References daughter_charge_, and mps_fire::i.

92 { return daughter_charge_[i]; }
std::vector< std::vector< int > > daughter_charge_
const std::vector<TLorentzVector>& reco::PFTau3ProngSummary::Daughter_P4 ( unsigned int  i) const
inline

Definition at line 93 of file PFTau3ProngSummary.h.

References daughter_p4_, and mps_fire::i.

93 { return daughter_p4_[i]; }
std::vector< std::vector< TLorentzVector > > daughter_p4_
const std::vector<int>& reco::PFTau3ProngSummary::Daughter_PDGID ( unsigned int  i) const
inline

Definition at line 91 of file PFTau3ProngSummary.h.

References daughter_PDGID_, and mps_fire::i.

91 { return daughter_PDGID_[i]; }
std::vector< std::vector< int > > daughter_PDGID_
CovMatrix reco::PFTau3ProngSummary::flightLenghtCov ( ) const
inline

Definition at line 70 of file PFTau3ProngSummary.h.

References TIP_.

70 { return TIP_->flightLengthCov(); }
reco::PFTauTransverseImpactParameterRef TIP_
const Vector& reco::PFTau3ProngSummary::flightLength ( ) const
inline

Definition at line 68 of file PFTau3ProngSummary.h.

References TIP_.

68 { return TIP_->flightLength(); }
reco::PFTauTransverseImpactParameterRef TIP_
double reco::PFTau3ProngSummary::flightLengthSig ( ) const
inline

Definition at line 69 of file PFTau3ProngSummary.h.

References TIP_.

69 { return TIP_->flightLengthSig(); }
reco::PFTauTransverseImpactParameterRef TIP_
bool reco::PFTau3ProngSummary::has3ProngSolution ( unsigned int  i) const
inline

Definition at line 84 of file PFTau3ProngSummary.h.

References has3ProngSolution_, and mps_fire::i.

Referenced by ~PFTau3ProngSummary().

84 { return has3ProngSolution_[i]; }
std::vector< bool > has3ProngSolution_
bool reco::PFTau3ProngSummary::hasSecondaryVertex ( ) const
inline

Definition at line 65 of file PFTau3ProngSummary.h.

References TIP_.

65 { return TIP_->hasSecondaryVertex(); }
reco::PFTauTransverseImpactParameterRef TIP_
const TVector3& reco::PFTau3ProngSummary::HelixFitSecondaryVertex ( ) const
inline

Definition at line 79 of file PFTau3ProngSummary.h.

References sv_.

79 { return sv_; }
const CovMatrix& reco::PFTau3ProngSummary::HelixFitSecondaryVertexCov ( ) const
inline

Definition at line 80 of file PFTau3ProngSummary.h.

References svcov_.

80 { return svcov_; }
double PFTau3ProngSummary::M_12 ( ) const

Definition at line 86 of file PFTau3ProngSummary.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, daughter_charge_, daughter_p4_, has3ProngSolution_, mps_fire::i, dqmiolumiharvest::j, and Tau_Charge().

Referenced by M_A1().

86  {
87  for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
88  if (has3ProngSolution_[i] == true) {
89  int charge = Tau_Charge();
90  TLorentzVector LV;
91  for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
92  if (daughter_charge_[i][j] == charge)
93  LV += daughter_p4_[i][j];
94  }
95  return LV.M();
96  }
97  }
98  return 0.0;
99 }
std::vector< bool > has3ProngSolution_
math::XYZTLorentzVectorD LV
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< TLorentzVector > > daughter_p4_
double PFTau3ProngSummary::M_13 ( ) const

Definition at line 100 of file PFTau3ProngSummary.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, daughter_charge_, daughter_p4_, newFWLiteAna::found, has3ProngSolution_, mps_fire::i, dqmiolumiharvest::j, and Tau_Charge().

Referenced by M_A1().

100  {
101  for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
102  if (has3ProngSolution_[i] == true) {
103  int charge = Tau_Charge();
104  TLorentzVector LV_opp;
105  for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
106  if (daughter_charge_[i][j] == -1 * charge)
107  LV_opp = daughter_p4_[i][j];
108  }
109  TLorentzVector LV_pair;
110  bool found(false);
111  for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
112  if (daughter_charge_[i][j] == charge) {
113  TLorentzVector LV = daughter_p4_[i][j];
114  LV += LV_opp;
115  if (!found)
116  LV_pair = LV;
117  else if (LV_pair.M() > LV.M())
118  LV_pair = LV;
119  found = true;
120  }
121  }
122  if (found)
123  return LV_pair.M();
124  }
125  }
126  return 0.0;
127 }
std::vector< bool > has3ProngSolution_
math::XYZTLorentzVectorD LV
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< TLorentzVector > > daughter_p4_
double PFTau3ProngSummary::M_23 ( ) const

Definition at line 129 of file PFTau3ProngSummary.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, daughter_charge_, daughter_p4_, newFWLiteAna::found, has3ProngSolution_, mps_fire::i, dqmiolumiharvest::j, and Tau_Charge().

Referenced by M_A1().

129  {
130  for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
131  if (has3ProngSolution_[i] == true) {
132  int charge = Tau_Charge();
133  TLorentzVector LV_opp;
134  for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
135  if (daughter_charge_[i][j] == -1 * charge)
136  LV_opp = daughter_p4_[i][j];
137  }
138  TLorentzVector LV_pair;
139  bool found(false);
140  for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
141  if (daughter_charge_[i][j] == charge) {
142  TLorentzVector LV = daughter_p4_[i][j];
143  LV += LV_opp;
144  if (!found)
145  LV_pair = LV;
146  else if (LV_pair.M() < LV.M())
147  LV_pair = LV;
148  found = true;
149  }
150  }
151  if (found)
152  return LV_pair.M();
153  }
154  }
155  return 0.0;
156 }
std::vector< bool > has3ProngSolution_
math::XYZTLorentzVectorD LV
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< TLorentzVector > > daughter_p4_
double reco::PFTau3ProngSummary::M_A1 ( ) const
inline

Definition at line 74 of file PFTau3ProngSummary.h.

References a1_, M_12(), M_13(), M_23(), and Tau_Charge().

74 { return a1_.M(); }
const reco::PFTauTransverseImpactParameterRef& reco::PFTau3ProngSummary::PFTauTIP ( ) const
inline

Definition at line 61 of file PFTau3ProngSummary.h.

References TIP_.

61 { return TIP_; }
reco::PFTauTransverseImpactParameterRef TIP_
const VertexRef& reco::PFTau3ProngSummary::primaryVertex ( ) const
inline

Definition at line 63 of file PFTau3ProngSummary.h.

References TIP_.

63 { return TIP_->primaryVertex(); }
reco::PFTauTransverseImpactParameterRef TIP_
CovMatrix reco::PFTau3ProngSummary::primaryVertexCov ( ) const
inline

Definition at line 64 of file PFTau3ProngSummary.h.

References TIP_.

64 { return TIP_->primaryVertexCov(); }
reco::PFTauTransverseImpactParameterRef TIP_
const VertexRef& reco::PFTau3ProngSummary::secondaryVertex ( ) const
inline

Definition at line 66 of file PFTau3ProngSummary.h.

References TIP_.

66 { return TIP_->secondaryVertex(); }
reco::PFTauTransverseImpactParameterRef TIP_
CovMatrix reco::PFTau3ProngSummary::secondaryVertexCov ( ) const
inline

Definition at line 67 of file PFTau3ProngSummary.h.

References TIP_.

67 { return TIP_->secondaryVertexCov(); }
reco::PFTauTransverseImpactParameterRef TIP_
double reco::PFTau3ProngSummary::SignificanceOfThetaGJ ( unsigned int  i) const
inline

Definition at line 86 of file PFTau3ProngSummary.h.

References mps_fire::i, and thetaGJsig_.

86  {
87  return thetaGJsig_[i];
88  } // 0 or less means the theta_GF has
std::vector< double > thetaGJsig_
double reco::PFTau3ProngSummary::Solution_Chi2 ( unsigned int  i) const
inline

Definition at line 85 of file PFTau3ProngSummary.h.

References mps_fire::i, and solution_Chi2_.

85 { return solution_Chi2_[i]; }
std::vector< double > solution_Chi2_
const TLorentzVector& reco::PFTau3ProngSummary::Tau ( unsigned int  i) const
inline

Definition at line 90 of file PFTau3ProngSummary.h.

References mps_fire::i, and tau_p4_.

90 { return tau_p4_[i]; }
std::vector< TLorentzVector > tau_p4_
int PFTau3ProngSummary::Tau_Charge ( ) const

Definition at line 158 of file PFTau3ProngSummary.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, daughter_charge_, daughter_p4_, has3ProngSolution_, mps_fire::i, and dqmiolumiharvest::j.

Referenced by M_12(), M_13(), M_23(), and M_A1().

158  {
159  for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
160  if (has3ProngSolution_[i] == true) {
161  int charge = 0;
162  for (unsigned int j = 0; j < daughter_p4_[i].size(); j++)
163  charge += daughter_charge_[i][j];
164  return charge;
165  }
166  }
167  return 0;
168 }
std::vector< bool > has3ProngSolution_
std::vector< std::vector< int > > daughter_charge_
std::vector< std::vector< TLorentzVector > > daughter_p4_
double reco::PFTau3ProngSummary::Vertex_chi2 ( ) const
inline

Definition at line 81 of file PFTau3ProngSummary.h.

References vertex_chi2_.

double reco::PFTau3ProngSummary::Vertex_ndf ( ) const
inline

Definition at line 82 of file PFTau3ProngSummary.h.

References vertex_ndf_.

82 { return vertex_ndf_; }
double reco::PFTau3ProngSummary::Vertex_Prob ( ) const
inline

Definition at line 83 of file PFTau3ProngSummary.h.

References vertex_chi2_, and vertex_ndf_.

Member Data Documentation

TLorentzVector reco::PFTau3ProngSummary::a1_
private

Definition at line 97 of file PFTau3ProngSummary.h.

Referenced by A1_LV(), M_A1(), and PFTau3ProngSummary().

std::vector<std::vector<int> > reco::PFTau3ProngSummary::daughter_charge_
private
std::vector<std::vector<TLorentzVector> > reco::PFTau3ProngSummary::daughter_p4_
private
std::vector<std::vector<int> > reco::PFTau3ProngSummary::daughter_PDGID_
private

Definition at line 106 of file PFTau3ProngSummary.h.

Referenced by AddSolution(), Daughter_PDGID(), and PFTau3ProngSummary().

std::vector<bool> reco::PFTau3ProngSummary::has3ProngSolution_
private
std::vector<double> reco::PFTau3ProngSummary::solution_Chi2_
private

Definition at line 103 of file PFTau3ProngSummary.h.

Referenced by AddSolution(), PFTau3ProngSummary(), and Solution_Chi2().

TVector3 reco::PFTau3ProngSummary::sv_
private

Definition at line 98 of file PFTau3ProngSummary.h.

Referenced by HelixFitSecondaryVertex(), and PFTau3ProngSummary().

CovMatrix reco::PFTau3ProngSummary::svcov_
private

Definition at line 99 of file PFTau3ProngSummary.h.

Referenced by HelixFitSecondaryVertexCov(), and PFTau3ProngSummary().

std::vector<TLorentzVector> reco::PFTau3ProngSummary::tau_p4_
private

Definition at line 105 of file PFTau3ProngSummary.h.

Referenced by AddSolution(), PFTau3ProngSummary(), and Tau().

std::vector<double> reco::PFTau3ProngSummary::thetaGJsig_
private

Definition at line 104 of file PFTau3ProngSummary.h.

Referenced by AddSolution(), PFTau3ProngSummary(), and SignificanceOfThetaGJ().

reco::PFTauTransverseImpactParameterRef reco::PFTau3ProngSummary::TIP_
private
double reco::PFTau3ProngSummary::vertex_chi2_
private

Definition at line 100 of file PFTau3ProngSummary.h.

Referenced by PFTau3ProngSummary(), Vertex_chi2(), and Vertex_Prob().

double reco::PFTau3ProngSummary::vertex_ndf_
private

Definition at line 101 of file PFTau3ProngSummary.h.

Referenced by PFTau3ProngSummary(), Vertex_ndf(), and Vertex_Prob().