CMS 3D CMS Logo

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