CMS 3D CMS Logo

TopDQMHelpers.cc
Go to the documentation of this file.
2 
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)
16  return massWBoson_;
17 }
18 
19 double Calculate::massTopQuark(const std::vector<reco::Jet>& jets) {
20  if (!failed_ && massTopQuark_ < 0)
22  return massTopQuark_;
23 }
24 
25 double Calculate::massBTopQuark(const std::vector<reco::Jet>& jets, std::vector<double> VbtagWP, double btagWP_) {
26  if (!failed_ && massBTopQuark_ < 0)
27  operator2(jets, VbtagWP, btagWP_);
28  return massBTopQuark_;
29 }
30 
32  if (tmassWBoson_ < 0)
33  operator()(b, mu, met);
34  return tmassWBoson_;
35 }
36 
38  if (tmassTopQuark_ < 0)
39  operator()(b, lepton, met);
40  return tmassTopQuark_;
41 }
42 
43 void Calculate::operator()(const reco::Jet& bJet, reco::RecoCandidate* lepton, const reco::MET& met) {
44  double metT = sqrt(pow(met.px(), 2) + pow(met.py(), 2));
45  double lepT = sqrt(pow(lepton->px(), 2) + pow(lepton->py(), 2));
46  double bT = sqrt(pow(bJet.px(), 2) + pow(bJet.py(), 2));
47  reco::Particle::LorentzVector WT = lepton->p4() + met.p4();
48  tmassWBoson_ = 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()) - (topT.py() * topT.py()));
51 }
52 
53 void Calculate::operator()(const std::vector<reco::Jet>& jets) {
54  if (maxNJets_ < 0)
55  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 = idx + 1; jdx < maxNJets_; ++jdx) {
71  //if (jdx <= idx) continue;
72  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
73  if (kdx == idx || kdx == jdx)
74  continue;
75  reco::Particle::LorentzVector sum = 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() + jets[maxPtIndices[2]].p4()).mass();
87 
88  // associate those jets that get closest to the W mass
89  // with their invariant mass to the W boson
90  double wDist = -1.;
91  std::vector<int> wMassIndices;
92  wMassIndices.push_back(-1);
93  wMassIndices.push_back(-1);
94  for (unsigned idx = 0; idx < maxPtIndices.size(); ++idx) {
95  for (unsigned jdx = 0; jdx < maxPtIndices.size(); ++jdx) {
96  if (jdx == idx || maxPtIndices[idx] > maxPtIndices[jdx])
97  continue;
98  reco::Particle::LorentzVector sum = jets[maxPtIndices[idx]].p4() + jets[maxPtIndices[jdx]].p4();
99  if (wDist < 0. || wDist > fabs(sum.mass() - wMass_)) {
100  wDist = fabs(sum.mass() - wMass_);
101  wMassIndices.clear();
102  wMassIndices.push_back(maxPtIndices[idx]);
103  wMassIndices.push_back(maxPtIndices[jdx]);
104  }
105  }
106  }
107  massWBoson_ = (jets[wMassIndices[0]].p4() + jets[wMassIndices[1]].p4()).mass();
108 }
109 
110 void Calculate::operator2(const std::vector<reco::Jet>& jets, std::vector<double> bjet, double btagWP) {
111  if (maxNJets_ < 0)
112  maxNJets_ = jets.size();
113  failed_ = jets.size() < (unsigned int)maxNJets_;
114  if (failed_) {
115  return;
116  }
117  if (jets.size() != bjet.size()) {
118  return;
119  }
120 
121  // associate those jets with maximum pt of the vectorial
122  // sum to the hadronic decay chain. Require ONLY 1 btagged jet
123  double maxBPt = -1.;
124  std::vector<int> maxBPtIndices;
125  maxBPtIndices.push_back(-1);
126  maxBPtIndices.push_back(-1);
127  maxBPtIndices.push_back(-1);
128  for (int idx = 0; idx < maxNJets_; ++idx) {
129  for (int jdx = idx + 1; jdx < maxNJets_; ++jdx) {
130  //if (jdx <= idx) continue;
131  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
132  if (kdx == idx || kdx == jdx)
133  continue;
134  // require only 1b-jet
135  if ((bjet[idx] > btagWP && bjet[jdx] <= btagWP && bjet[kdx] <= btagWP) ||
136  (bjet[idx] <= btagWP && bjet[jdx] > btagWP && bjet[kdx] <= btagWP) ||
137  (bjet[idx] <= btagWP && bjet[jdx] <= btagWP && bjet[kdx] > btagWP)) {
138  reco::Particle::LorentzVector sum = jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
139  if (maxBPt < 0. || maxBPt < sum.pt()) {
140  maxBPt = sum.pt();
141  maxBPtIndices.clear();
142  maxBPtIndices.push_back(idx);
143  maxBPtIndices.push_back(jdx);
144  maxBPtIndices.push_back(kdx);
145  }
146  }
147  }
148  }
149  }
150  if (maxBPtIndices[0] < 0 || maxBPtIndices[1] < 0 || maxBPtIndices[2] < 0)
151  return;
152  massBTopQuark_ = (jets[maxBPtIndices[0]].p4() + jets[maxBPtIndices[1]].p4() + jets[maxBPtIndices[2]].p4()).mass();
153 }
154 
156  : failed_(false),
157  maxNJets_(maxNJets),
158  wMass_(wMass),
159  massWBoson_(-1.),
160  massTopQuark_(-1.),
161  massBTopQuark_(-1.),
162  tmassWBoson_(-1),
163  tmassTopQuark_(-1) {}
164 
165 double Calculate_miniAOD::massWBoson(const std::vector<pat::Jet>& jets) {
166  if (!failed_ && massWBoson_ < 0)
167  operator()(jets);
168  return massWBoson_;
169 }
170 
171 double Calculate_miniAOD::massTopQuark(const std::vector<pat::Jet>& jets) {
172  if (!failed_ && massTopQuark_ < 0)
173  operator()(jets);
174  return massTopQuark_;
175 }
176 
177 double Calculate_miniAOD::massBTopQuark(const std::vector<pat::Jet>& jets,
178  std::vector<double> VbtagWP,
179  double btagWP_) {
180  if (!failed_ && massBTopQuark_ < 0)
181  operator2(jets, VbtagWP, btagWP_);
182  return massBTopQuark_;
183 }
184 
186  if (tmassWBoson_ < 0)
187  operator()(b, mu, met);
188  return tmassWBoson_;
189 }
190 
192  if (tmassWBoson_ < 0)
193  operator()(b, mu, met);
194  return tmassWBoson_;
195 }
196 
198  if (tmassTopQuark_ < 0)
199  operator()(b, lepton, met);
200  return tmassTopQuark_;
201 }
202 
204  if (tmassTopQuark_ < 0)
205  operator()(b, lepton, met);
206  return tmassTopQuark_;
207 }
208 
209 void Calculate_miniAOD::operator()(const pat::Jet& bJet, pat::Muon* lepton, const pat::MET& met) {
210  double metT = sqrt(pow(met.px(), 2) + pow(met.py(), 2));
211  double lepT = sqrt(pow(lepton->px(), 2) + pow(lepton->py(), 2));
212  double bT = sqrt(pow(bJet.px(), 2) + pow(bJet.py(), 2));
213  reco::Particle::LorentzVector WT = lepton->p4() + met.p4();
214  tmassWBoson_ = sqrt(pow(metT + lepT, 2) - (WT.px() * WT.px()) - (WT.py() * WT.py()));
215  reco::Particle::LorentzVector topT = WT + bJet.p4();
216  tmassTopQuark_ = sqrt(pow((metT + lepT + bT), 2) - (topT.px() * topT.px()) - (topT.py() * topT.py()));
217 }
218 
220  double metT = sqrt(pow(met.px(), 2) + pow(met.py(), 2));
221  double lepT = sqrt(pow(lepton->px(), 2) + pow(lepton->py(), 2));
222  double bT = sqrt(pow(bJet.px(), 2) + pow(bJet.py(), 2));
223  reco::Particle::LorentzVector WT = lepton->p4() + met.p4();
224  tmassWBoson_ = sqrt(pow(metT + lepT, 2) - (WT.px() * WT.px()) - (WT.py() * WT.py()));
225  reco::Particle::LorentzVector topT = WT + bJet.p4();
226  tmassTopQuark_ = sqrt(pow((metT + lepT + bT), 2) - (topT.px() * topT.px()) - (topT.py() * topT.py()));
227 }
228 
229 void Calculate_miniAOD::operator()(const std::vector<pat::Jet>& jets) {
230  if (maxNJets_ < 0)
231  maxNJets_ = jets.size();
232  failed_ = jets.size() < (unsigned int)maxNJets_;
233  if (failed_) {
234  return;
235  }
236 
237  // associate those jets with maximum pt of the vectorial
238  // sum to the hadronic decay chain
239  double maxPt = -1.;
240  std::vector<int> maxPtIndices;
241  maxPtIndices.push_back(-1);
242  maxPtIndices.push_back(-1);
243  maxPtIndices.push_back(-1);
244 
245  for (int idx = 0; idx < maxNJets_; ++idx) {
246  for (int jdx = idx + 1; jdx < maxNJets_; ++jdx) {
247  //if (jdx <= idx) continue;
248  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
249  if (kdx == idx || kdx == jdx)
250  continue;
251  reco::Particle::LorentzVector sum = jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
252  if (maxPt < 0. || maxPt < sum.pt()) {
253  maxPt = sum.pt();
254  maxPtIndices.clear();
255  maxPtIndices.push_back(idx);
256  maxPtIndices.push_back(jdx);
257  maxPtIndices.push_back(kdx);
258  }
259  }
260  }
261  }
262  massTopQuark_ = (jets[maxPtIndices[0]].p4() + jets[maxPtIndices[1]].p4() + jets[maxPtIndices[2]].p4()).mass();
263 
264  // associate those jets that get closest to the W mass
265  // with their invariant mass to the W boson
266  double wDist = -1.;
267  std::vector<int> wMassIndices;
268  wMassIndices.push_back(-1);
269  wMassIndices.push_back(-1);
270  for (unsigned idx = 0; idx < maxPtIndices.size(); ++idx) {
271  for (unsigned jdx = 0; jdx < maxPtIndices.size(); ++jdx) {
272  if (jdx == idx || maxPtIndices[idx] > maxPtIndices[jdx])
273  continue;
274  reco::Particle::LorentzVector sum = jets[maxPtIndices[idx]].p4() + jets[maxPtIndices[jdx]].p4();
275  if (wDist < 0. || wDist > fabs(sum.mass() - wMass_)) {
276  wDist = fabs(sum.mass() - wMass_);
277  wMassIndices.clear();
278  wMassIndices.push_back(maxPtIndices[idx]);
279  wMassIndices.push_back(maxPtIndices[jdx]);
280  }
281  }
282  }
283  massWBoson_ = (jets[wMassIndices[0]].p4() + jets[wMassIndices[1]].p4()).mass();
284 }
285 
286 void Calculate_miniAOD::operator2(const std::vector<pat::Jet>& jets, std::vector<double> bjet, double btagWP) {
287  if (maxNJets_ < 0)
288  maxNJets_ = jets.size();
289  failed_ = jets.size() < (unsigned int)maxNJets_;
290 
291  if (failed_) {
292  return;
293  }
294 
295  if (jets.size() != bjet.size()) {
296  return;
297  }
298 
299  // associate those jets with maximum pt of the vectorial
300  // sum to the hadronic decay chain. Require ONLY 1 btagged jet
301  double maxBPt = -1.;
302  std::vector<int> maxBPtIndices;
303  maxBPtIndices.push_back(-1);
304  maxBPtIndices.push_back(-1);
305  maxBPtIndices.push_back(-1);
306  for (int idx = 0; idx < maxNJets_; ++idx) {
307  for (int jdx = idx + 1; jdx < maxNJets_; ++jdx) {
308  //if (jdx <= idx) continue;
309  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
310  if (kdx == idx || kdx == jdx)
311  continue;
312  // require only 1b-jet
313  if ((bjet[idx] > btagWP && bjet[jdx] <= btagWP && bjet[kdx] <= btagWP) ||
314  (bjet[idx] <= btagWP && bjet[jdx] > btagWP && bjet[kdx] <= btagWP) ||
315  (bjet[idx] <= btagWP && bjet[jdx] <= btagWP && bjet[kdx] > btagWP)) {
316  reco::Particle::LorentzVector sum = jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
317  if (maxBPt < 0. || maxBPt < sum.pt()) {
318  maxBPt = sum.pt();
319  maxBPtIndices.clear();
320  maxBPtIndices.push_back(idx);
321  maxBPtIndices.push_back(jdx);
322  maxBPtIndices.push_back(kdx);
323  }
324  }
325  }
326  }
327  }
328  if (maxBPtIndices[0] < 0 || maxBPtIndices[1] < 0 || maxBPtIndices[2] < 0)
329  return;
330  massBTopQuark_ = (jets[maxBPtIndices[0]].p4() + jets[maxBPtIndices[1]].p4() + jets[maxBPtIndices[2]].p4()).mass();
331 }
Calculate_miniAOD::Calculate_miniAOD
Calculate_miniAOD(int maxNJets, double wMass)
default constructor
Definition: TopDQMHelpers.cc:155
Calculate_miniAOD::tmassTopQuark
double tmassTopQuark(pat::Muon *lep, const pat::MET &met, const pat::Jet &b)
calculate top quark transverse mass estimate
Definition: TopDQMHelpers.cc:203
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
funct::false
false
Definition: Factorize.h:29
Calculate::massBTopQuark_
double massBTopQuark_
cache of b-tagged top quark mass estimate
Definition: TopDQMHelpers.h:170
amptDefaultParameters_cff.mu
mu
Definition: amptDefaultParameters_cff.py:16
Calculate::Calculate
Calculate(int maxNJets, double wMass)
default constructor
Definition: TopDQMHelpers.cc:3
Calculate_miniAOD::maxNJets_
int maxNJets_
max. number of jets to be considered
Definition: TopDQMHelpers.h:112
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
pat::Muon
Analysis-level muon class.
Definition: Muon.h:51
Calculate::failed_
bool failed_
indicate failed associations
Definition: TopDQMHelpers.h:160
reco::Particle::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
BTaggingMonitor_cfi.met
met
Definition: BTaggingMonitor_cfi.py:84
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
MuonErrorMatrixAnalyzer_cfi.maxPt
maxPt
Definition: MuonErrorMatrixAnalyzer_cfi.py:19
reco::MET
Definition: MET.h:41
pat::Jet
Analysis-level calorimeter jet class.
Definition: Jet.h:77
reco::LeafCandidate::py
double py() const final
y coordinate of momentum vector
Definition: LeafCandidate.h:142
Calculate_miniAOD::operator()
void operator()(const std::vector< pat::Jet > &jets)
Definition: TopDQMHelpers.cc:229
Calculate::operator()
void operator()(const std::vector< reco::Jet > &jets)
Definition: TopDQMHelpers.cc:53
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
Calculate_miniAOD::wMass_
double wMass_
paramater of the w boson mass
Definition: TopDQMHelpers.h:114
Calculate_miniAOD::operator2
void operator2(const std::vector< pat::Jet > &, std::vector< double >, double)
do the calculation of the t-quark mass with one b-jet
Definition: TopDQMHelpers.cc:286
Calculate::operator2
void operator2(const std::vector< reco::Jet > &, std::vector< double >, double)
do the calculation of the t-quark mass with one b-jet
Definition: TopDQMHelpers.cc:110
Calculate::massTopQuark_
double massTopQuark_
cache of top quark mass estimate
Definition: TopDQMHelpers.h:168
b
double b
Definition: hdecay.h:118
Calculate::maxNJets_
int maxNJets_
max. number of jets to be considered
Definition: TopDQMHelpers.h:162
Calculate_miniAOD::massWBoson_
double massWBoson_
cache of w boson mass estimate
Definition: TopDQMHelpers.h:116
Calculate_miniAOD::massWBoson
double massWBoson(const std::vector< pat::Jet > &jets)
calculate W boson mass estimate
Definition: TopDQMHelpers.cc:165
createfilelist.int
int
Definition: createfilelist.py:10
reco::LeafCandidate::p4
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:114
SUSYDQMAnalyzer_cfi.maxNJets
maxNJets
Definition: SUSYDQMAnalyzer_cfi.py:13
pat::MET
Analysis-level MET class.
Definition: MET.h:40
Calculate::massWBoson_
double massWBoson_
cache of w boson mass estimate
Definition: TopDQMHelpers.h:166
Calculate::tmassWBoson_
double tmassWBoson_
cache of W boson transverse mass estimate
Definition: TopDQMHelpers.h:172
Calculate_miniAOD::tmassTopQuark_
double tmassTopQuark_
cache of top quark transverse mass estimate
Definition: TopDQMHelpers.h:124
Calculate_miniAOD::massBTopQuark_
double massBTopQuark_
cache of b-tagged top quark mass estimate
Definition: TopDQMHelpers.h:120
reco::RecoCandidate
Definition: RecoCandidate.h:20
Calculate::tmassTopQuark_
double tmassTopQuark_
cache of top quark transverse mass estimate
Definition: TopDQMHelpers.h:174
reco::GsfElectron::p4
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:211
Calculate_miniAOD::massTopQuark
double massTopQuark(const std::vector< pat::Jet > &jets)
calculate t-quark mass estimate
Definition: TopDQMHelpers.cc:171
Calculate::massTopQuark
double massTopQuark(const std::vector< reco::Jet > &jets)
calculate t-quark mass estimate
Definition: TopDQMHelpers.cc:19
Calculate::tmassWBoson
double tmassWBoson(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate W boson transverse mass estimate
Definition: TopDQMHelpers.cc:31
Calculate::tmassTopQuark
double tmassTopQuark(reco::RecoCandidate *lep, const reco::MET &met, const reco::Jet &b)
calculate top quark transverse mass estimate
Definition: TopDQMHelpers.cc:37
Calculate_miniAOD::tmassWBoson_
double tmassWBoson_
cache of W boson transverse mass estimate
Definition: TopDQMHelpers.h:122
Calculate::massBTopQuark
double massBTopQuark(const std::vector< reco::Jet > &jets, std::vector< double > VbtagWP, double btagWP_)
calculate b-tagged t-quark mass estimate
Definition: TopDQMHelpers.cc:25
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
pseudoTop_cfi.wMass
wMass
Definition: pseudoTop_cfi.py:12
Calculate::massWBoson
double massWBoson(const std::vector< reco::Jet > &jets)
calculate W boson mass estimate
Definition: TopDQMHelpers.cc:13
Calculate::wMass_
double wMass_
paramater of the w boson mass
Definition: TopDQMHelpers.h:164
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Calculate_miniAOD::massBTopQuark
double massBTopQuark(const std::vector< pat::Jet > &jets, std::vector< double > VbtagWP, double btagWP_)
calculate b-tagged t-quark mass estimate
Definition: TopDQMHelpers.cc:177
pat::Electron
Analysis-level electron class.
Definition: Electron.h:51
reco::LeafCandidate::px
double px() const final
x coordinate of momentum vector
Definition: LeafCandidate.h:140
Calculate_miniAOD::failed_
bool failed_
indicate failed associations
Definition: TopDQMHelpers.h:110
Calculate_miniAOD::massTopQuark_
double massTopQuark_
cache of top quark mass estimate
Definition: TopDQMHelpers.h:118
TopDQMHelpers.h
Calculate_miniAOD::tmassWBoson
double tmassWBoson(pat::Muon *lep, const pat::MET &met, const pat::Jet &b)
calculate W boson transverse mass estimate
Definition: TopDQMHelpers.cc:185