File indexing completed on 2024-04-06 12:05:17
0001 #include "DataFormats/TauReco/interface/PFTau3ProngSummary.h"
0002 #include "TMatrixT.h"
0003 #include "TMatrixTSym.h"
0004 #include "TVectorT.h"
0005 using namespace reco;
0006
0007 PFTau3ProngSummary::PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,
0008 TLorentzVector a1,
0009 double vertex_chi2,
0010 double vertex_ndf) {
0011 TIP_ = TIP;
0012 for (unsigned int i = 0; i < nsolutions; i++) {
0013 has3ProngSolution_.push_back(false);
0014 solution_Chi2_.push_back(0);
0015 thetaGJsig_.push_back(0);
0016 tau_p4_.push_back(TLorentzVector(0, 0, 0, 0));
0017 daughter_PDGID_.push_back(std::vector<int>());
0018 daughter_charge_.push_back(std::vector<int>());
0019 daughter_p4_.push_back(std::vector<TLorentzVector>());
0020 }
0021 a1_ = a1;
0022 sv_ = TVector3(TIP_->secondaryVertex()->x(), TIP_->secondaryVertex()->y(), TIP_->secondaryVertex()->z());
0023 svcov_ = TIP_->secondaryVertexCov();
0024 vertex_chi2_ = vertex_chi2;
0025 vertex_ndf_ = vertex_ndf;
0026 }
0027
0028 PFTau3ProngSummary::PFTau3ProngSummary(reco::PFTauTransverseImpactParameterRef TIP,
0029 TLorentzVector a1,
0030 double vertex_chi2,
0031 double vertex_ndf,
0032 TVector3 sv,
0033 CovMatrix svcov) {
0034 TIP_ = TIP;
0035 for (unsigned int i = 0; i < nsolutions; i++) {
0036 has3ProngSolution_.push_back(false);
0037 solution_Chi2_.push_back(0);
0038 thetaGJsig_.push_back(0);
0039 tau_p4_.push_back(TLorentzVector(0, 0, 0, 0));
0040 daughter_PDGID_.push_back(std::vector<int>());
0041 daughter_charge_.push_back(std::vector<int>());
0042 daughter_p4_.push_back(std::vector<TLorentzVector>());
0043 }
0044 a1_ = a1;
0045 sv_ = sv;
0046 svcov_ = svcov;
0047 vertex_chi2_ = vertex_chi2;
0048 vertex_ndf_ = vertex_ndf;
0049 }
0050
0051 PFTau3ProngSummary::PFTau3ProngSummary() {
0052 for (unsigned int i = 0; i < nsolutions; i++) {
0053 has3ProngSolution_.push_back(false);
0054 solution_Chi2_.push_back(0);
0055 thetaGJsig_.push_back(0);
0056 tau_p4_.push_back(TLorentzVector(0, 0, 0, 0));
0057 daughter_PDGID_.push_back(std::vector<int>());
0058 daughter_charge_.push_back(std::vector<int>());
0059 daughter_p4_.push_back(std::vector<TLorentzVector>());
0060 }
0061 }
0062
0063 PFTau3ProngSummary* PFTau3ProngSummary::clone() const { return new PFTau3ProngSummary(*this); }
0064
0065 bool PFTau3ProngSummary::AddSolution(unsigned int solution,
0066 const TLorentzVector& tau,
0067 const std::vector<TLorentzVector>& daughter_p4,
0068 const std::vector<int>& daughter_charge,
0069 const std::vector<int>& daughter_PDGID,
0070 bool has3ProngSolution,
0071 double solutionChi2,
0072 double thetaGJsig) {
0073 if (solution < nsolutions) {
0074 has3ProngSolution_[solution] = true;
0075 solution_Chi2_[solution] = solutionChi2;
0076 thetaGJsig_[solution] = thetaGJsig;
0077 tau_p4_[solution] = tau;
0078 daughter_PDGID_[solution] = daughter_PDGID;
0079 daughter_charge_[solution] = daughter_charge;
0080 daughter_p4_[solution] = daughter_p4;
0081 return true;
0082 }
0083 return false;
0084 }
0085
0086 double PFTau3ProngSummary::M_12() const {
0087 for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
0088 if (has3ProngSolution_[i] == true) {
0089 int charge = Tau_Charge();
0090 TLorentzVector LV;
0091 for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
0092 if (daughter_charge_[i][j] == charge)
0093 LV += daughter_p4_[i][j];
0094 }
0095 return LV.M();
0096 }
0097 }
0098 return 0.0;
0099 }
0100 double PFTau3ProngSummary::M_13() const {
0101 for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
0102 if (has3ProngSolution_[i] == true) {
0103 int charge = Tau_Charge();
0104 TLorentzVector LV_opp;
0105 for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
0106 if (daughter_charge_[i][j] == -1 * charge)
0107 LV_opp = daughter_p4_[i][j];
0108 }
0109 TLorentzVector LV_pair;
0110 bool found(false);
0111 for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
0112 if (daughter_charge_[i][j] == charge) {
0113 TLorentzVector LV = daughter_p4_[i][j];
0114 LV += LV_opp;
0115 if (!found)
0116 LV_pair = LV;
0117 else if (LV_pair.M() > LV.M())
0118 LV_pair = LV;
0119 found = true;
0120 }
0121 }
0122 if (found)
0123 return LV_pair.M();
0124 }
0125 }
0126 return 0.0;
0127 }
0128
0129 double PFTau3ProngSummary::M_23() const {
0130 for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
0131 if (has3ProngSolution_[i] == true) {
0132 int charge = Tau_Charge();
0133 TLorentzVector LV_opp;
0134 for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
0135 if (daughter_charge_[i][j] == -1 * charge)
0136 LV_opp = daughter_p4_[i][j];
0137 }
0138 TLorentzVector LV_pair;
0139 bool found(false);
0140 for (unsigned int j = 0; j < daughter_p4_[i].size(); j++) {
0141 if (daughter_charge_[i][j] == charge) {
0142 TLorentzVector LV = daughter_p4_[i][j];
0143 LV += LV_opp;
0144 if (!found)
0145 LV_pair = LV;
0146 else if (LV_pair.M() < LV.M())
0147 LV_pair = LV;
0148 found = true;
0149 }
0150 }
0151 if (found)
0152 return LV_pair.M();
0153 }
0154 }
0155 return 0.0;
0156 }
0157
0158 int PFTau3ProngSummary::Tau_Charge() const {
0159 for (unsigned int i = 0; i < has3ProngSolution_.size(); i++) {
0160 if (has3ProngSolution_[i] == true) {
0161 int charge = 0;
0162 for (unsigned int j = 0; j < daughter_p4_[i].size(); j++)
0163 charge += daughter_charge_[i][j];
0164 return charge;
0165 }
0166 }
0167 return 0;
0168 }