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 Calculate::Calculate(int maxNJets, double wMass):
4  failed_(false), maxNJets_(maxNJets), wMass_(wMass), massWBoson_(-1.), massTopQuark_(-1.)
5 {
6 }
7 
8 double
9 Calculate::massWBoson(const std::vector<reco::Jet>& jets)
10 {
11  if(!failed_&& massWBoson_<0) operator()(jets); return massWBoson_;
12 }
13 
14 double
15 Calculate::massTopQuark(const std::vector<reco::Jet>& jets)
16 {
17  if(!failed_&& massTopQuark_<0) operator()(jets); return massTopQuark_;
18 }
19 
20 void
21 Calculate::operator()(const std::vector<reco::Jet>& jets)
22 {
23  if(maxNJets_<0) maxNJets_=jets.size();
24  failed_= jets.size()<(unsigned int) maxNJets_;
25  if( failed_){ return; }
26 
27  // associate those jets with maximum pt of the vectorial
28  // sum to the hadronic decay chain
29  double maxPt=-1.;
30  std::vector<int> maxPtIndices;
31  maxPtIndices.push_back(-1);
32  maxPtIndices.push_back(-1);
33  maxPtIndices.push_back(-1);
34  for(int idx=0; idx<maxNJets_; ++idx){
35  for(int jdx=0; jdx<maxNJets_; ++jdx){ if(jdx<=idx) continue;
36  for(int kdx=0; kdx<maxNJets_; ++kdx){ if(kdx==idx || kdx==jdx) continue;
37  reco::Particle::LorentzVector sum = jets[idx].p4()+jets[jdx].p4()+jets[kdx].p4();
38  if( maxPt<0. || maxPt<sum.pt() ){
39  maxPt=sum.pt();
40  maxPtIndices.clear();
41  maxPtIndices.push_back(idx);
42  maxPtIndices.push_back(jdx);
43  maxPtIndices.push_back(kdx);
44  }
45  }
46  }
47  }
48  massTopQuark_= (jets[maxPtIndices[0]].p4()+
49  jets[maxPtIndices[1]].p4()+
50  jets[maxPtIndices[2]].p4()).mass();
51 
52  // associate those jets that get closest to the W mass
53  // with their invariant mass to the W boson
54  double wDist =-1.;
55  std::vector<int> wMassIndices;
56  wMassIndices.push_back(-1);
57  wMassIndices.push_back(-1);
58  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){
59  for(unsigned jdx=0; jdx<maxPtIndices.size(); ++jdx){
60  if( jdx==idx || maxPtIndices[idx]>maxPtIndices[jdx] ) continue;
61  reco::Particle::LorentzVector sum = jets[maxPtIndices[idx]].p4()+jets[maxPtIndices[jdx]].p4();
62  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
63  wDist=fabs(sum.mass()-wMass_);
64  wMassIndices.clear();
65  wMassIndices.push_back(maxPtIndices[idx]);
66  wMassIndices.push_back(maxPtIndices[jdx]);
67  }
68  }
69  }
70  massWBoson_= (jets[wMassIndices[0]].p4()+
71  jets[wMassIndices[1]].p4()).mass();
72 }
double massWBoson_
cache of w boson mass estimate
Definition: TopDQMHelpers.h:91
double massTopQuark_
cache of top quark mass estimate
Definition: TopDQMHelpers.h:93
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
Definition: TopDQMHelpers.cc:9
int maxNJets_
max. number of jets to be considered
Definition: TopDQMHelpers.h:87
vector< PseudoJet > jets
double wMass_
paramater of the w boson mass
Definition: TopDQMHelpers.h:89
void operator()(const std::vector< reco::Jet > &jets)
bool failed_
indicate failed associations
Definition: TopDQMHelpers.h:85
tuple mass
Definition: scaleCards.py:27
Calculate(int maxNJets, double wMass)
default constructor
Definition: TopDQMHelpers.cc:3
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25