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),
5  maxNJets_(maxNJets),
6  wMass_(wMass),
7  massWBoson_(-1.),
8  massTopQuark_(-1.),
9  massBTopQuark_(-1.),
10  tmassWBoson_(-1),
11  tmassTopQuark_(-1) {}
12 
13 double Calculate::massWBoson(const std::vector<reco::Jet>& jets) {
14  if (!failed_ && massWBoson_ < 0) operator()(jets);
15  return massWBoson_;
16 }
17 
18 double Calculate::massTopQuark(const std::vector<reco::Jet>& jets) {
19  if (!failed_ && massTopQuark_ < 0) operator()(jets);
20  return massTopQuark_;
21 }
22 
23 double Calculate::massBTopQuark(const std::vector<reco::Jet>& jets,
24  std::vector<double> VbtagWP, double btagWP_) {
25  if (!failed_ && massBTopQuark_ < 0) operator2(jets, VbtagWP, btagWP_);
26  return massBTopQuark_;
27 }
28 
30  const reco::Jet& b) {
31  if (tmassWBoson_ < 0) operator()(b, mu, met);
32  return tmassWBoson_;
33 }
34 
36  const reco::MET& met, const reco::Jet& b) {
37  if (tmassTopQuark_ < 0) operator()(b, lepton, met);
38  return tmassTopQuark_;
39 }
40 
42  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_ =
48  sqrt(pow(metT + lepT, 2) - (WT.px() * WT.px()) - (WT.py() * WT.py()));
49  reco::Particle::LorentzVector topT = WT + bJet.p4();
50  tmassTopQuark_ = sqrt(pow((metT + lepT + bT), 2) - (topT.px() * topT.px()) -
51  (topT.py() * topT.py()));
52 }
53 
54 void Calculate::operator()(const std::vector<reco::Jet>& jets) {
55  if (maxNJets_ < 0) maxNJets_ = jets.size();
56  failed_ = jets.size() < (unsigned int)maxNJets_;
57  if (failed_) {
58  return;
59  }
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) {
71  if (jdx <= idx) continue;
72  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
73  if (kdx == idx || kdx == jdx) continue;
75  jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
76  if (maxPt < 0. || maxPt < sum.pt()) {
77  maxPt = sum.pt();
78  maxPtIndices.clear();
79  maxPtIndices.push_back(idx);
80  maxPtIndices.push_back(jdx);
81  maxPtIndices.push_back(kdx);
82  }
83  }
84  }
85  }
86  massTopQuark_ = (jets[maxPtIndices[0]].p4() + jets[maxPtIndices[1]].p4() +
87  jets[maxPtIndices[2]].p4()).mass();
88 
89  // associate those jets that get closest to the W mass
90  // with their invariant mass to the W boson
91  double wDist = -1.;
92  std::vector<int> wMassIndices;
93  wMassIndices.push_back(-1);
94  wMassIndices.push_back(-1);
95  for (unsigned idx = 0; idx < maxPtIndices.size(); ++idx) {
96  for (unsigned jdx = 0; jdx < maxPtIndices.size(); ++jdx) {
97  if (jdx == idx || maxPtIndices[idx] > maxPtIndices[jdx]) continue;
99  jets[maxPtIndices[idx]].p4() + jets[maxPtIndices[jdx]].p4();
100  if (wDist < 0. || wDist > fabs(sum.mass() - wMass_)) {
101  wDist = fabs(sum.mass() - wMass_);
102  wMassIndices.clear();
103  wMassIndices.push_back(maxPtIndices[idx]);
104  wMassIndices.push_back(maxPtIndices[jdx]);
105  }
106  }
107  }
108  massWBoson_ =
109  (jets[wMassIndices[0]].p4() + jets[wMassIndices[1]].p4()).mass();
110 }
111 
112 void Calculate::operator2(const std::vector<reco::Jet>& jets,
113  std::vector<double> bjet, double btagWP) {
114  if (maxNJets_ < 0) maxNJets_ = jets.size();
115  failed_ = jets.size() < (unsigned int)maxNJets_;
116  if (failed_) {
117  return;
118  }
119  if (jets.size() != bjet.size()) {
120  return;
121  }
122 
123  // associate those jets with maximum pt of the vectorial
124  // sum to the hadronic decay chain. Require ONLY 1 btagged jet
125  double maxBPt = -1.;
126  std::vector<int> maxBPtIndices;
127  maxBPtIndices.push_back(-1);
128  maxBPtIndices.push_back(-1);
129  maxBPtIndices.push_back(-1);
130  for (int idx = 0; idx < maxNJets_; ++idx) {
131  for (int jdx = 0; jdx < maxNJets_; ++jdx) {
132  if (jdx <= idx) continue;
133  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
134  if (kdx == idx || kdx == jdx) continue;
135  // require only 1b-jet
136  if ((bjet[idx] > btagWP && bjet[jdx] <= btagWP &&
137  bjet[kdx] <= btagWP) ||
138  (bjet[idx] <= btagWP && bjet[jdx] > btagWP &&
139  bjet[kdx] <= btagWP) ||
140  (bjet[idx] <= btagWP && bjet[jdx] <= btagWP &&
141  bjet[kdx] > btagWP)) {
143  jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
144  if (maxBPt < 0. || maxBPt < sum.pt()) {
145  maxBPt = sum.pt();
146  maxBPtIndices.clear();
147  maxBPtIndices.push_back(idx);
148  maxBPtIndices.push_back(jdx);
149  maxBPtIndices.push_back(kdx);
150  }
151  }
152  }
153  }
154  }
155  if (maxBPtIndices[0] < 0 || maxBPtIndices[1] < 0 || maxBPtIndices[2] < 0)
156  return;
157  massBTopQuark_ = (jets[maxBPtIndices[0]].p4() + jets[maxBPtIndices[1]].p4() +
158  jets[maxBPtIndices[2]].p4()).mass();
159 }
double massWBoson_
cache of w boson mass estimate
double massBTopQuark(const std::vector< reco::Jet > &jets, std::vector< double > VbtagWP, double btagWP_)
calculate b-tagged t-quark mass estimate
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
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:42
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
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
virtual double px() const
x coordinate of momentum vector
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:3
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
virtual double py() const
y coordinate of momentum vector
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