CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
Calculate_miniAOD Class Reference

#include <TopDQMHelpers.h>

Public Member Functions

 Calculate_miniAOD (int maxNJets, double wMass)
 default constructor More...
 
double massBTopQuark (const std::vector< pat::Jet > &jets, std::vector< double > VbtagWP, double btagWP_)
 calculate b-tagged t-quark mass estimate More...
 
double massTopQuark (const std::vector< pat::Jet > &jets)
 calculate t-quark mass estimate More...
 
double massWBoson (const std::vector< pat::Jet > &jets)
 calculate W boson mass estimate More...
 
double tmassTopQuark (pat::Muon *lep, const pat::MET &met, const pat::Jet &b)
 calculate top quark transverse mass estimate More...
 
double tmassTopQuark (pat::Electron *lep, const pat::MET &met, const pat::Jet &b)
 calculate top quark transverse mass estimate More...
 
double tmassWBoson (pat::Muon *lep, const pat::MET &met, const pat::Jet &b)
 calculate W boson transverse mass estimate More...
 
double tmassWBoson (pat::Electron *lep, const pat::MET &met, const pat::Jet &b)
 calculate W boson transverse mass estimate More...
 
 ~Calculate_miniAOD ()
 default destructor More...
 

Private Member Functions

void operator() (const std::vector< pat::Jet > &jets)
 
void operator() (const pat::Jet &bJet, pat::Electron *lepton, const pat::MET &met)
 do the calculation of the transverse top and W masses More...
 
void operator() (const pat::Jet &bJet, pat::Muon *lepton, const pat::MET &met)
 
void operator2 (const std::vector< pat::Jet > &, std::vector< double >, double)
 do the calculation of the t-quark mass with one b-jet More...
 

Private Attributes

bool failed_
 indicate failed associations More...
 
double massBTopQuark_
 cache of b-tagged top quark mass estimate More...
 
double massTopQuark_
 cache of top quark mass estimate More...
 
double massWBoson_
 cache of w boson mass estimate More...
 
int maxNJets_
 max. number of jets to be considered More...
 
double tmassTopQuark_
 cache of top quark transverse mass estimate More...
 
double tmassWBoson_
 cache of W boson transverse mass estimate More...
 
double wMass_
 paramater of the w boson mass More...
 

Detailed Description

Definition at line 72 of file TopDQMHelpers.h.

Constructor & Destructor Documentation

Calculate_miniAOD::Calculate_miniAOD ( int  maxNJets,
double  wMass 
)

default constructor

Definition at line 168 of file TopDQMHelpers.cc.

169  : failed_(false),
171  wMass_(wMass),
172  massWBoson_(-1.),
173  massTopQuark_(-1.),
174  massBTopQuark_(-1.),
175  tmassWBoson_(-1),
176  tmassTopQuark_(-1) {}
double massBTopQuark_
cache of b-tagged top quark mass estimate
double massTopQuark_
cache of top quark mass estimate
bool failed_
indicate failed associations
double massWBoson_
cache of w boson mass estimate
double tmassWBoson_
cache of W boson transverse mass estimate
int maxNJets_
max. number of jets to be considered
double wMass_
paramater of the w boson mass
double tmassTopQuark_
cache of top quark transverse mass estimate
Calculate_miniAOD::~Calculate_miniAOD ( )
inline

default destructor

Definition at line 77 of file TopDQMHelpers.h.

77 {};

Member Function Documentation

double Calculate_miniAOD::massBTopQuark ( const std::vector< pat::Jet > &  jets,
std::vector< double >  VbtagWP,
double  btagWP_ 
)

calculate b-tagged t-quark mass estimate

Definition at line 188 of file TopDQMHelpers.cc.

References failed_, massBTopQuark_, and operator2().

Referenced by TopSingleLepton_miniAOD::MonitorEnsemble::fill().

189  {
190 
191  if (!failed_ && massBTopQuark_ < 0) operator2(jets, VbtagWP, btagWP_);
192  return massBTopQuark_;
193 }
double massBTopQuark_
cache of b-tagged top quark mass estimate
bool failed_
indicate failed associations
vector< PseudoJet > jets
void operator2(const std::vector< pat::Jet > &, std::vector< double >, double)
do the calculation of the t-quark mass with one b-jet
double Calculate_miniAOD::massTopQuark ( const std::vector< pat::Jet > &  jets)

calculate t-quark mass estimate

Definition at line 183 of file TopDQMHelpers.cc.

References failed_, massTopQuark_, and operator()().

Referenced by SingleTopTChannelLepton_miniAOD::MonitorEnsemble::fill(), and TopSingleLepton_miniAOD::MonitorEnsemble::fill().

183  {
184  if (!failed_ && massTopQuark_ < 0) operator()(jets);
185  return massTopQuark_;
186 }
double massTopQuark_
cache of top quark mass estimate
bool failed_
indicate failed associations
void operator()(const std::vector< pat::Jet > &jets)
vector< PseudoJet > jets
double Calculate_miniAOD::massWBoson ( const std::vector< pat::Jet > &  jets)

calculate W boson mass estimate

Definition at line 178 of file TopDQMHelpers.cc.

References failed_, massWBoson_, and operator()().

Referenced by SingleTopTChannelLepton_miniAOD::MonitorEnsemble::fill(), and TopSingleLepton_miniAOD::MonitorEnsemble::fill().

178  {
179  if (!failed_ && massWBoson_ < 0) operator()(jets);
180  return massWBoson_;
181 }
bool failed_
indicate failed associations
double massWBoson_
cache of w boson mass estimate
void operator()(const std::vector< pat::Jet > &jets)
vector< PseudoJet > jets
void Calculate_miniAOD::operator() ( const std::vector< pat::Jet > &  jets)
private

do the calculation; this is called only once per event by the first function call to return a mass estimate. The once calculated values are cached afterwards

Definition at line 250 of file TopDQMHelpers.cc.

References failed_, customizeTrackingMonitorSeedNumber::idx, massTopQuark_, massWBoson_, maxNJets_, and wMass_.

Referenced by massTopQuark(), massWBoson(), tmassTopQuark(), and tmassWBoson().

250  {
251  if (maxNJets_ < 0) maxNJets_ = jets.size();
252  failed_ = jets.size() < (unsigned int)maxNJets_;
253  if (failed_) {
254  return;
255  }
256 
257  // associate those jets with maximum pt of the vectorial
258  // sum to the hadronic decay chain
259  double maxPt = -1.;
260  std::vector<int> maxPtIndices;
261  maxPtIndices.push_back(-1);
262  maxPtIndices.push_back(-1);
263  maxPtIndices.push_back(-1);
264 
265  for (int idx = 0; idx < maxNJets_; ++idx) {
266  for (int jdx = idx+1; jdx < maxNJets_; ++jdx) {
267  //if (jdx <= idx) continue;
268  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
269  if (kdx == idx || kdx == jdx) continue;
271  jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
272  if (maxPt < 0. || maxPt < sum.pt()) {
273  maxPt = sum.pt();
274  maxPtIndices.clear();
275  maxPtIndices.push_back(idx);
276  maxPtIndices.push_back(jdx);
277  maxPtIndices.push_back(kdx);
278  }
279  }
280  }
281  }
282  massTopQuark_ = (jets[maxPtIndices[0]].p4() + jets[maxPtIndices[1]].p4() +
283  jets[maxPtIndices[2]].p4()).mass();
284 
285  // associate those jets that get closest to the W mass
286  // with their invariant mass to the W boson
287  double wDist = -1.;
288  std::vector<int> wMassIndices;
289  wMassIndices.push_back(-1);
290  wMassIndices.push_back(-1);
291  for (unsigned idx = 0; idx < maxPtIndices.size(); ++idx) {
292  for (unsigned jdx = 0; jdx < maxPtIndices.size(); ++jdx) {
293  if (jdx == idx || maxPtIndices[idx] > maxPtIndices[jdx]) continue;
295  jets[maxPtIndices[idx]].p4() + jets[maxPtIndices[jdx]].p4();
296  if (wDist < 0. || wDist > fabs(sum.mass() - wMass_)) {
297  wDist = fabs(sum.mass() - wMass_);
298  wMassIndices.clear();
299  wMassIndices.push_back(maxPtIndices[idx]);
300  wMassIndices.push_back(maxPtIndices[jdx]);
301  }
302  }
303  }
304  massWBoson_ =
305  (jets[wMassIndices[0]].p4() + jets[wMassIndices[1]].p4()).mass();
306 }
double massTopQuark_
cache of top quark mass estimate
bool failed_
indicate failed associations
double massWBoson_
cache of w boson mass estimate
vector< PseudoJet > jets
int maxNJets_
max. number of jets to be considered
double wMass_
paramater of the w boson mass
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
void Calculate_miniAOD::operator() ( const pat::Jet bJet,
pat::Electron lepton,
const pat::MET met 
)
private

do the calculation of the transverse top and W masses

Definition at line 236 of file TopDQMHelpers.cc.

References reco::LeafCandidate::p4(), reco::GsfElectron::p4(), funct::pow(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), mathSSE::sqrt(), tmassTopQuark_, and tmassWBoson_.

237  {
238  double metT = sqrt(pow(met.px(), 2) + pow(met.py(), 2));
239  double lepT = sqrt(pow(lepton->px(), 2) + pow(lepton->py(), 2));
240  double bT = sqrt(pow(bJet.px(), 2) + pow(bJet.py(), 2));
241  reco::Particle::LorentzVector WT = lepton->p4() + met.p4();
242  tmassWBoson_ =
243  sqrt(pow(metT + lepT, 2) - (WT.px() * WT.px()) - (WT.py() * WT.py()));
244  reco::Particle::LorentzVector topT = WT + bJet.p4();
245  tmassTopQuark_ = sqrt(pow((metT + lepT + bT), 2) - (topT.px() * topT.px()) -
246  (topT.py() * topT.py()));
247 }
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:223
double tmassWBoson_
cache of W boson transverse mass estimate
T sqrt(T t)
Definition: SSEVec.h:18
virtual double py() const final
y coordinate of momentum vector
double tmassTopQuark_
cache of top quark transverse mass estimate
virtual double px() const final
x coordinate of momentum vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void Calculate_miniAOD::operator() ( const pat::Jet bJet,
pat::Muon lepton,
const pat::MET met 
)
private

Definition at line 222 of file TopDQMHelpers.cc.

References reco::LeafCandidate::p4(), funct::pow(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), mathSSE::sqrt(), tmassTopQuark_, and tmassWBoson_.

223  {
224  double metT = sqrt(pow(met.px(), 2) + pow(met.py(), 2));
225  double lepT = sqrt(pow(lepton->px(), 2) + pow(lepton->py(), 2));
226  double bT = sqrt(pow(bJet.px(), 2) + pow(bJet.py(), 2));
227  reco::Particle::LorentzVector WT = lepton->p4() + met.p4();
228  tmassWBoson_ =
229  sqrt(pow(metT + lepT, 2) - (WT.px() * WT.px()) - (WT.py() * WT.py()));
230  reco::Particle::LorentzVector topT = WT + bJet.p4();
231  tmassTopQuark_ = sqrt(pow((metT + lepT + bT), 2) - (topT.px() * topT.px()) -
232  (topT.py() * topT.py()));
233 }
double tmassWBoson_
cache of W boson transverse mass estimate
T sqrt(T t)
Definition: SSEVec.h:18
virtual double py() const final
y coordinate of momentum vector
double tmassTopQuark_
cache of top quark transverse mass estimate
virtual double px() const final
x coordinate of momentum vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void Calculate_miniAOD::operator2 ( const std::vector< pat::Jet > &  jets,
std::vector< double >  bjet,
double  btagWP 
)
private

do the calculation of the t-quark mass with one b-jet

Definition at line 308 of file TopDQMHelpers.cc.

References failed_, customizeTrackingMonitorSeedNumber::idx, massBTopQuark_, and maxNJets_.

Referenced by massBTopQuark().

309  {
310  if (maxNJets_ < 0) maxNJets_ = jets.size();
311  failed_ = jets.size() < (unsigned int)maxNJets_;
312 
313  if (failed_) {
314  return;
315  }
316 
317  if (jets.size() != bjet.size()) {
318  return;
319  }
320 
321  // associate those jets with maximum pt of the vectorial
322  // sum to the hadronic decay chain. Require ONLY 1 btagged jet
323  double maxBPt = -1.;
324  std::vector<int> maxBPtIndices;
325  maxBPtIndices.push_back(-1);
326  maxBPtIndices.push_back(-1);
327  maxBPtIndices.push_back(-1);
328  for (int idx = 0; idx < maxNJets_; ++idx) {
329  for (int jdx = idx+1; jdx < maxNJets_; ++jdx) {
330  //if (jdx <= idx) continue;
331  for (int kdx = 0; kdx < maxNJets_; ++kdx) {
332  if (kdx == idx || kdx == jdx) continue;
333  // require only 1b-jet
334  if ((bjet[idx] > btagWP && bjet[jdx] <= btagWP &&
335  bjet[kdx] <= btagWP) ||
336  (bjet[idx] <= btagWP && bjet[jdx] > btagWP &&
337  bjet[kdx] <= btagWP) ||
338  (bjet[idx] <= btagWP && bjet[jdx] <= btagWP &&
339  bjet[kdx] > btagWP)) {
341  jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
342  if (maxBPt < 0. || maxBPt < sum.pt()) {
343  maxBPt = sum.pt();
344  maxBPtIndices.clear();
345  maxBPtIndices.push_back(idx);
346  maxBPtIndices.push_back(jdx);
347  maxBPtIndices.push_back(kdx);
348  }
349  }
350  }
351  }
352  }
353  if (maxBPtIndices[0] < 0 || maxBPtIndices[1] < 0 || maxBPtIndices[2] < 0)
354  return;
355  massBTopQuark_ = (jets[maxBPtIndices[0]].p4() + jets[maxBPtIndices[1]].p4() +
356  jets[maxBPtIndices[2]].p4()).mass();
357 }
double massBTopQuark_
cache of b-tagged top quark mass estimate
bool failed_
indicate failed associations
vector< PseudoJet > jets
int maxNJets_
max. number of jets to be considered
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
double Calculate_miniAOD::tmassTopQuark ( pat::Muon lep,
const pat::MET met,
const pat::Jet b 
)

calculate top quark transverse mass estimate

Definition at line 215 of file TopDQMHelpers.cc.

References operator()(), and tmassTopQuark_.

Referenced by SingleTopTChannelLepton_miniAOD::MonitorEnsemble::fill().

216  {
217  if (tmassTopQuark_ < 0) operator()(b, lepton, met);
218  return tmassTopQuark_;
219 }
void operator()(const std::vector< pat::Jet > &jets)
double tmassTopQuark_
cache of top quark transverse mass estimate
double Calculate_miniAOD::tmassTopQuark ( pat::Electron lep,
const pat::MET met,
const pat::Jet b 
)

calculate top quark transverse mass estimate

Definition at line 208 of file TopDQMHelpers.cc.

References operator()(), and tmassTopQuark_.

209  {
210  if (tmassTopQuark_ < 0) operator()(b, lepton, met);
211  return tmassTopQuark_;
212 }
void operator()(const std::vector< pat::Jet > &jets)
double tmassTopQuark_
cache of top quark transverse mass estimate
double Calculate_miniAOD::tmassWBoson ( pat::Muon lep,
const pat::MET met,
const pat::Jet b 
)

calculate W boson transverse mass estimate

Definition at line 195 of file TopDQMHelpers.cc.

References operator()(), and tmassWBoson_.

Referenced by SingleTopTChannelLepton_miniAOD::MonitorEnsemble::fill().

196  {
197  if (tmassWBoson_ < 0) operator()(b, mu, met);
198  return tmassWBoson_;
199 }
double tmassWBoson_
cache of W boson transverse mass estimate
void operator()(const std::vector< pat::Jet > &jets)
const int mu
Definition: Constants.h:22
double Calculate_miniAOD::tmassWBoson ( pat::Electron lep,
const pat::MET met,
const pat::Jet b 
)

calculate W boson transverse mass estimate

Definition at line 201 of file TopDQMHelpers.cc.

References operator()(), and tmassWBoson_.

202  {
203  if (tmassWBoson_ < 0) operator()(b, mu, met);
204  return tmassWBoson_;
205 }
double tmassWBoson_
cache of W boson transverse mass estimate
void operator()(const std::vector< pat::Jet > &jets)
const int mu
Definition: Constants.h:22

Member Data Documentation

bool Calculate_miniAOD::failed_
private

indicate failed associations

Definition at line 120 of file TopDQMHelpers.h.

Referenced by massBTopQuark(), massTopQuark(), massWBoson(), operator()(), and operator2().

double Calculate_miniAOD::massBTopQuark_
private

cache of b-tagged top quark mass estimate

Definition at line 130 of file TopDQMHelpers.h.

Referenced by massBTopQuark(), and operator2().

double Calculate_miniAOD::massTopQuark_
private

cache of top quark mass estimate

Definition at line 128 of file TopDQMHelpers.h.

Referenced by massTopQuark(), and operator()().

double Calculate_miniAOD::massWBoson_
private

cache of w boson mass estimate

Definition at line 126 of file TopDQMHelpers.h.

Referenced by massWBoson(), and operator()().

int Calculate_miniAOD::maxNJets_
private

max. number of jets to be considered

Definition at line 122 of file TopDQMHelpers.h.

Referenced by operator()(), and operator2().

double Calculate_miniAOD::tmassTopQuark_
private

cache of top quark transverse mass estimate

Definition at line 134 of file TopDQMHelpers.h.

Referenced by operator()(), and tmassTopQuark().

double Calculate_miniAOD::tmassWBoson_
private

cache of W boson transverse mass estimate

Definition at line 132 of file TopDQMHelpers.h.

Referenced by operator()(), and tmassWBoson().

double Calculate_miniAOD::wMass_
private

paramater of the w boson mass

Definition at line 124 of file TopDQMHelpers.h.

Referenced by operator()().