CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFTau3ProngSummary.cc
Go to the documentation of this file.
2 #include "TMatrixT.h"
3 #include "TMatrixTSym.h"
5 #include "TVectorT.h"
6 using namespace reco;
7 
8 PFTau3ProngSummary::PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,TLorentzVector a1,double vertex_chi2,double vertex_ndf){
9  TIP_=TIP;
10  for(unsigned int i=0;i<nsolutions;i++){
11  has3ProngSolution_.push_back(false);
12  solution_Chi2_.push_back(0);
13  thetaGJsig_.push_back(0);
14  tau_p4_.push_back(TLorentzVector(0,0,0,0));
15  daughter_PDGID_.push_back(std::vector<int>());
16  daughter_charge_.push_back(std::vector<int>());
17  daughter_p4_.push_back(std::vector<TLorentzVector>());
18  }
19  a1_=a1;
20  sv_=TVector3(TIP_->secondaryVertex()->x(),TIP_->secondaryVertex()->y(),TIP_->secondaryVertex()->z());
21  svcov_=TIP_->secondaryVertexCov();
22  vertex_chi2_=vertex_chi2;
23  vertex_ndf_=vertex_ndf;
24 }
25 
26 PFTau3ProngSummary::PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,TLorentzVector a1,double vertex_chi2,double vertex_ndf,TVector3 sv,CovMatrix svcov){
27  TIP_=TIP;
28  for(unsigned int i=0;i<nsolutions;i++){
29  has3ProngSolution_.push_back(false);
30  solution_Chi2_.push_back(0);
31  thetaGJsig_.push_back(0);
32  tau_p4_.push_back(TLorentzVector(0,0,0,0));
33  daughter_PDGID_.push_back(std::vector<int>());
34  daughter_charge_.push_back(std::vector<int>());
35  daughter_p4_.push_back(std::vector<TLorentzVector>());
36  }
37  a1_=a1;
38  sv_=sv;
39  svcov_=svcov;
40  vertex_chi2_=vertex_chi2;
41  vertex_ndf_=vertex_ndf;
42 }
43 
44 
46  for(unsigned int i=0;i<nsolutions;i++){
47  has3ProngSolution_.push_back(false);
48  solution_Chi2_.push_back(0);
49  thetaGJsig_.push_back(0);
50  tau_p4_.push_back(TLorentzVector(0,0,0,0));
51  daughter_PDGID_.push_back(std::vector<int>());
52  daughter_charge_.push_back(std::vector<int>());
53  daughter_p4_.push_back(std::vector<TLorentzVector>());
54  }
55 }
56 
58  return new PFTau3ProngSummary(*this);
59 }
60 
61 
62 bool PFTau3ProngSummary::AddSolution(unsigned int solution, const TLorentzVector& tau, const std::vector<TLorentzVector>& daughter_p4,
63  const std::vector<int>& daughter_charge, const std::vector<int>& daughter_PDGID,
64  bool has3ProngSolution, double solutionChi2, double thetaGJsig){
65  if(solution<nsolutions){
66  has3ProngSolution_[solution]=true;
67  solution_Chi2_[solution]=solutionChi2;
68  thetaGJsig_[solution]=thetaGJsig;
69  tau_p4_[solution]=tau;
70  daughter_PDGID_[solution]=daughter_PDGID;
71  daughter_charge_[solution]=daughter_charge;
72  daughter_p4_[solution]=daughter_p4;
73  return true;
74  }
75  return false;
76 }
77 
78 
80  for(unsigned int i=0;i<has3ProngSolution_.size();i++){
81  if(has3ProngSolution_[i]==true){
82  int charge=Tau_Charge();
83  TLorentzVector LV;
84  for(unsigned int j=0;j<daughter_p4_[i].size();j++){
85  if(daughter_charge_[i][j]==charge)LV+=daughter_p4_[i][j];
86  }
87  return LV.M();
88  }
89  }
90  return 0.0;
91 }
93  for(unsigned int i=0;i<has3ProngSolution_.size();i++){
94  if(has3ProngSolution_[i]==true){
95  int charge=Tau_Charge();
96  TLorentzVector LV_opp;
97  for(unsigned int j=0;j<daughter_p4_[i].size();j++){
98  if(daughter_charge_[i][j]==-1*charge)LV_opp=daughter_p4_[i][j];
99  }
100  TLorentzVector LV_pair;
101  bool found(false);
102  for(unsigned int j=0;j<daughter_p4_[i].size();j++){
103  if(daughter_charge_[i][j]==charge){
104  TLorentzVector LV=daughter_p4_[i][j];
105  LV+=LV_opp;
106  if(!found)LV_pair=LV;
107  else if(LV_pair.M()>LV.M())LV_pair=LV;
108  found=true;
109  }
110  }
111  if(found)return LV_pair.M();
112  }
113  }
114  return 0.0;
115 }
116 
118  for(unsigned int i=0;i<has3ProngSolution_.size();i++){
119  if(has3ProngSolution_[i]==true){
120  int charge=Tau_Charge();
121  TLorentzVector LV_opp;
122  for(unsigned int j=0;j<daughter_p4_[i].size();j++){
123  if(daughter_charge_[i][j]==-1*charge)LV_opp=daughter_p4_[i][j];
124  }
125  TLorentzVector LV_pair;
126  bool found(false);
127  for(unsigned int j=0;j<daughter_p4_[i].size();j++){
128  if(daughter_charge_[i][j]==charge){
129  TLorentzVector LV=daughter_p4_[i][j];
130  LV+=LV_opp;
131  if(!found)LV_pair=LV;
132  else if(LV_pair.M()<LV.M())LV_pair=LV;
133  found=true;
134  }
135  }
136  if(found)return LV_pair.M();
137  }
138  }
139  return 0.0;
140 }
141 
143  for(unsigned int i=0;i<has3ProngSolution_.size();i++){
144  if(has3ProngSolution_[i]==true){
145  int charge;
146  for(unsigned int j=0;j<daughter_p4_[i].size();j++)charge+=daughter_charge_[i][j];
147  return charge;
148  }
149  }
150  return 0;
151 }
int i
Definition: DBlmapReader.cc:9
reco::PFTauTransverseImpactParameterRef TIP_
math::Error< dimension >::type CovMatrix
std::vector< bool > has3ProngSolution_
std::vector< double > thetaGJsig_
std::vector< double > solution_Chi2_
double charge(const std::vector< uint8_t > &Ampls)
math::XYZTLorentzVectorD LV
int j
Definition: DBlmapReader.cc:9
std::vector< TLorentzVector > tau_p4_
std::vector< std::vector< int > > daughter_charge_
PFTau3ProngSummary()
constructor from values
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)
PFTau3ProngSummary * clone() const
std::vector< std::vector< int > > daughter_PDGID_
std::vector< std::vector< TLorentzVector > > daughter_p4_