CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TopDQMHelpers.cc
Go to the documentation of this file.
2 
3 
4 Calculate::Calculate(int maxNJets, double wMass):
5  failed_(false), maxNJets_(maxNJets), wMass_(wMass), massWBoson_(-1.), massTopQuark_(-1.), massBTopQuark_(-1.), tmassWBoson_(-1),tmassTopQuark_(-1)
6 {
7 }
8 
9 double
10 Calculate::massWBoson(const std::vector<reco::Jet>& jets)
11 {
12  if(!failed_&& massWBoson_<0) operator()(jets); return massWBoson_;
13 }
14 
15 double
16 Calculate::massTopQuark(const std::vector<reco::Jet>& jets)
17 {
18  if(!failed_&& massTopQuark_<0) operator()(jets); return massTopQuark_;
19 }
20 
21 
22 double
23 Calculate::massBTopQuark(const std::vector<reco::Jet>& jets, std::vector<double> VbtagWP, double btagWP_)
24 {
25  if(!failed_&& massBTopQuark_<0) operator2(jets, VbtagWP, btagWP_); return massBTopQuark_;
26 }
27 
28 double
30 {
31  if( tmassWBoson_<0) operator()(b,mu,met); return tmassWBoson_;
32 }
33 
34 
35 double
37 {
38  if( tmassTopQuark_<0) operator()(b,lepton,met); return tmassTopQuark_;
39 }
40 
41 
42 void Calculate::operator()( const reco::Jet& bJet, reco::RecoCandidate* lepton, const reco::MET& met){
43  double metT = sqrt(pow(met.px(),2) + pow(met.py(),2));
44  double lepT = sqrt(pow(lepton->px(),2) + pow(lepton->py(),2));
45  double bT = sqrt(pow(bJet.px(),2) + pow(bJet.py(),2));
46  reco::Particle::LorentzVector WT = lepton->p4() + met.p4();
47  tmassWBoson_ = sqrt(pow(metT+lepT,2) - (WT.px()*WT.px()) - (WT.py()*WT.py()));
48  reco::Particle::LorentzVector topT = WT + bJet.p4();
49  tmassTopQuark_ = sqrt(pow((metT+lepT+bT),2) - (topT.px()*topT.px()) - (topT.py()*topT.py()));
50 }
51 
52 
53 
54 void
55 Calculate::operator()(const std::vector<reco::Jet>& jets)
56 {
57  if(maxNJets_<0) maxNJets_=jets.size();
58  failed_= jets.size()<(unsigned int) maxNJets_;
59  if( failed_){ return; }
60 
61  // associate those jets with maximum pt of the vectorial
62  // sum to the hadronic decay chain
63  double maxPt=-1.;
64  std::vector<int> maxPtIndices;
65  maxPtIndices.push_back(-1);
66  maxPtIndices.push_back(-1);
67  maxPtIndices.push_back(-1);
68 
69  for(int idx=0; idx<maxNJets_; ++idx){
70  for(int jdx=0; jdx<maxNJets_; ++jdx){ if(jdx<=idx) continue;
71  for(int kdx=0; kdx<maxNJets_; ++kdx){ if(kdx==idx || kdx==jdx) continue;
72  reco::Particle::LorentzVector sum = jets[idx].p4()+jets[jdx].p4()+jets[kdx].p4();
73  if( maxPt<0. || maxPt<sum.pt() ){
74  maxPt=sum.pt();
75  maxPtIndices.clear();
76  maxPtIndices.push_back(idx);
77  maxPtIndices.push_back(jdx);
78  maxPtIndices.push_back(kdx);
79  }
80  }
81  }
82  }
83  massTopQuark_= (jets[maxPtIndices[0]].p4()+
84  jets[maxPtIndices[1]].p4()+
85  jets[maxPtIndices[2]].p4()).mass();
86 
87  // associate those jets that get closest to the W mass
88  // with their invariant mass to the W boson
89  double wDist =-1.;
90  std::vector<int> wMassIndices;
91  wMassIndices.push_back(-1);
92  wMassIndices.push_back(-1);
93  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){
94  for(unsigned jdx=0; jdx<maxPtIndices.size(); ++jdx){
95  if( jdx==idx || maxPtIndices[idx]>maxPtIndices[jdx] ) continue;
96  reco::Particle::LorentzVector sum = jets[maxPtIndices[idx]].p4()+jets[maxPtIndices[jdx]].p4();
97  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
98  wDist=fabs(sum.mass()-wMass_);
99  wMassIndices.clear();
100  wMassIndices.push_back(maxPtIndices[idx]);
101  wMassIndices.push_back(maxPtIndices[jdx]);
102  }
103  }
104  }
105  massWBoson_= (jets[wMassIndices[0]].p4()+
106  jets[wMassIndices[1]].p4()).mass();
107 }
108 
109 void
110 Calculate::operator2(const std::vector<reco::Jet>& jets, std::vector<double> bjet, double btagWP)
111 {
112  if(maxNJets_<0) maxNJets_=jets.size();
113  failed_= jets.size()<(unsigned int) maxNJets_;
114  if( failed_){ return; }
115  if (jets.size() != bjet.size()){return;}
116 
117  // associate those jets with maximum pt of the vectorial
118  // sum to the hadronic decay chain. Require ONLY 1 btagged jet
119  double maxBPt=-1.;
120  std::vector<int> maxBPtIndices;
121  maxBPtIndices.push_back(-1);
122  maxBPtIndices.push_back(-1);
123  maxBPtIndices.push_back(-1);
124  for(int idx=0; idx<maxNJets_; ++idx){
125  for(int jdx=0; jdx<maxNJets_; ++jdx){ if(jdx<=idx) continue;
126  for(int kdx=0; kdx<maxNJets_; ++kdx){if(kdx==idx || kdx==jdx) continue;
127  //require only 1b-jet
128  if ((bjet[idx]> btagWP && bjet[jdx]<= btagWP && bjet[kdx]<= btagWP) ||
129  (bjet[idx]<= btagWP && bjet[jdx]> btagWP && bjet[kdx]<= btagWP) ||
130  (bjet[idx]<= btagWP && bjet[jdx]<= btagWP && bjet[kdx]> btagWP) ){
131  reco::Particle::LorentzVector sum = jets[idx].p4()+jets[jdx].p4()+jets[kdx].p4();
132  if( maxBPt<0. || maxBPt<sum.pt() ){
133  maxBPt=sum.pt();
134  maxBPtIndices.clear();
135  maxBPtIndices.push_back(idx);
136  maxBPtIndices.push_back(jdx);
137  maxBPtIndices.push_back(kdx);
138  }
139  }
140  }
141  }
142  }
143  if (maxBPtIndices[0]<0 || maxBPtIndices[1]<0 || maxBPtIndices[2]<0) return;
144  massBTopQuark_= (jets[maxBPtIndices[0]].p4()+
145  jets[maxBPtIndices[1]].p4()+
146  jets[maxBPtIndices[2]].p4()).mass();
147 }
double massWBoson_
cache of w boson mass estimate
tuple met
____________________________________________________________________________||
Definition: CaloMET_cfi.py:7
double massBTopQuark(const std::vector< reco::Jet > &jets, std::vector< double > VbtagWP, double btagWP_)
calculate b-tagged t-quark mass estimate
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
double massTopQuark_
cache of top quark mass estimate
double tmassTopQuark_
cache of top quark transverse mass estimate
Base class for all types of Jets.
Definition: Jet.h:20
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate t-quark mass estimate
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
void operator2(const std::vector< reco::Jet > &, std::vector< double >, double)
do the calculation of the t-quark mass with one b-jet
int maxNJets_
max. number of jets to be considered
Definition: MET.h:32
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
const int mu
Definition: Constants.h:22
double wMass_
paramater of the w boson mass
double tmassTopQuark(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate top quark transverse mass estimate
double tmassWBoson(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate W boson transverse mass estimate
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double b
Definition: hdecay.h:120
void operator()(const std::vector< reco::Jet > &jets)
bool failed_
indicate failed associations
volatile std::atomic< bool > shutdown_flag false
Calculate(int maxNJets, double wMass)
default constructor
Definition: TopDQMHelpers.cc:4
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double tmassWBoson_
cache of W boson transverse mass estimate
double massBTopQuark_
cache of b-tagged top quark mass estimate