CMS 3D CMS Logo

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 71 of file TopDQMHelpers.h.

Constructor & Destructor Documentation

◆ Calculate_miniAOD()

Calculate_miniAOD::Calculate_miniAOD ( int  maxNJets,
double  wMass 
)

default constructor

Definition at line 155 of file TopDQMHelpers.cc.

156  : failed_(false),
158  wMass_(wMass),
159  massWBoson_(-1.),
160  massTopQuark_(-1.),
161  massBTopQuark_(-1.),
162  tmassWBoson_(-1),
163  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::~Calculate_miniAOD ( )
inline

default destructor

Definition at line 76 of file TopDQMHelpers.h.

76 {}

Member Function Documentation

◆ massBTopQuark()

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 177 of file TopDQMHelpers.cc.

References failed_, PDWG_EXODelayedJetMET_cff::jets, massBTopQuark_, and operator2().

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

179  {
180  if (!failed_ && massBTopQuark_ < 0)
181  operator2(jets, VbtagWP, btagWP_);
182  return massBTopQuark_;
183 }
double massBTopQuark_
cache of b-tagged top quark mass estimate
bool failed_
indicate failed associations
void operator2(const std::vector< pat::Jet > &, std::vector< double >, double)
do the calculation of the t-quark mass with one b-jet

◆ massTopQuark()

double Calculate_miniAOD::massTopQuark ( const std::vector< pat::Jet > &  jets)

calculate t-quark mass estimate

Definition at line 171 of file TopDQMHelpers.cc.

References failed_, PDWG_EXODelayedJetMET_cff::jets, massTopQuark_, and operator()().

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

171  {
172  if (!failed_ && massTopQuark_ < 0)
173  operator()(jets);
174  return massTopQuark_;
175 }
double massTopQuark_
cache of top quark mass estimate
bool failed_
indicate failed associations
void operator()(const std::vector< pat::Jet > &jets)

◆ massWBoson()

double Calculate_miniAOD::massWBoson ( const std::vector< pat::Jet > &  jets)

calculate W boson mass estimate

Definition at line 165 of file TopDQMHelpers.cc.

References failed_, PDWG_EXODelayedJetMET_cff::jets, massWBoson_, and operator()().

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

165  {
166  if (!failed_ && massWBoson_ < 0)
167  operator()(jets);
168  return massWBoson_;
169 }
bool failed_
indicate failed associations
double massWBoson_
cache of w boson mass estimate
void operator()(const std::vector< pat::Jet > &jets)

◆ operator()() [1/3]

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 229 of file TopDQMHelpers.cc.

References failed_, heavyIonCSV_trainingSettings::idx, createfilelist::int, PDWG_EXODelayedJetMET_cff::jets, EgHLTOffHistBins_cfi::mass, massTopQuark_, massWBoson_, maxNJets_, PV_cfg::maxPt, and wMass_.

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

229  {
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 }
double massTopQuark_
cache of top quark mass estimate
bool failed_
indicate failed associations
double massWBoson_
cache of w boson mass estimate
int maxNJets_
max. number of jets to be considered
double wMass_
paramater of the w boson mass
maxPt
Definition: PV_cfg.py:224
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21

◆ operator()() [2/3]

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 219 of file TopDQMHelpers.cc.

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

219  {
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 }
const LorentzVector & p4() const final
four-momentum Lorentz vector
double px() const final
x coordinate of momentum vector
double tmassWBoson_
cache of W boson transverse mass estimate
T sqrt(T t)
Definition: SSEVec.h:23
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:217
double py() const final
y coordinate of momentum vector
double tmassTopQuark_
cache of top quark transverse mass estimate
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ operator()() [3/3]

void Calculate_miniAOD::operator() ( const pat::Jet bJet,
pat::Muon lepton,
const pat::MET met 
)
private

Definition at line 209 of file TopDQMHelpers.cc.

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

209  {
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 }
const LorentzVector & p4() const final
four-momentum Lorentz vector
double px() const final
x coordinate of momentum vector
double tmassWBoson_
cache of W boson transverse mass estimate
T sqrt(T t)
Definition: SSEVec.h:23
double py() const final
y coordinate of momentum vector
double tmassTopQuark_
cache of top quark transverse mass estimate
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ operator2()

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 286 of file TopDQMHelpers.cc.

References failed_, heavyIonCSV_trainingSettings::idx, createfilelist::int, PDWG_EXODelayedJetMET_cff::jets, EgHLTOffHistBins_cfi::mass, massBTopQuark_, and maxNJets_.

Referenced by massBTopQuark().

286  {
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 }
double massBTopQuark_
cache of b-tagged top quark mass estimate
bool failed_
indicate failed associations
int maxNJets_
max. number of jets to be considered
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21

◆ tmassTopQuark() [1/2]

double Calculate_miniAOD::tmassTopQuark ( pat::Muon lep,
const pat::MET met,
const pat::Jet b 
)

calculate top quark transverse mass estimate

Definition at line 203 of file TopDQMHelpers.cc.

References b, BTaggingMonitor_cfi::met, operator()(), and tmassTopQuark_.

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

203  {
204  if (tmassTopQuark_ < 0)
205  operator()(b, lepton, met);
206  return tmassTopQuark_;
207 }
void operator()(const std::vector< pat::Jet > &jets)
double b
Definition: hdecay.h:120
double tmassTopQuark_
cache of top quark transverse mass estimate

◆ tmassTopQuark() [2/2]

double Calculate_miniAOD::tmassTopQuark ( pat::Electron lep,
const pat::MET met,
const pat::Jet b 
)

calculate top quark transverse mass estimate

Definition at line 197 of file TopDQMHelpers.cc.

References b, BTaggingMonitor_cfi::met, operator()(), and tmassTopQuark_.

197  {
198  if (tmassTopQuark_ < 0)
199  operator()(b, lepton, met);
200  return tmassTopQuark_;
201 }
void operator()(const std::vector< pat::Jet > &jets)
double b
Definition: hdecay.h:120
double tmassTopQuark_
cache of top quark transverse mass estimate

◆ tmassWBoson() [1/2]

double Calculate_miniAOD::tmassWBoson ( pat::Muon lep,
const pat::MET met,
const pat::Jet b 
)

calculate W boson transverse mass estimate

Definition at line 185 of file TopDQMHelpers.cc.

References b, BTaggingMonitor_cfi::met, amptDefaultParameters_cff::mu, operator()(), and tmassWBoson_.

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

185  {
186  if (tmassWBoson_ < 0)
187  operator()(b, mu, met);
188  return tmassWBoson_;
189 }
double tmassWBoson_
cache of W boson transverse mass estimate
void operator()(const std::vector< pat::Jet > &jets)
double b
Definition: hdecay.h:120

◆ tmassWBoson() [2/2]

double Calculate_miniAOD::tmassWBoson ( pat::Electron lep,
const pat::MET met,
const pat::Jet b 
)

calculate W boson transverse mass estimate

Definition at line 191 of file TopDQMHelpers.cc.

References b, BTaggingMonitor_cfi::met, amptDefaultParameters_cff::mu, operator()(), and tmassWBoson_.

191  {
192  if (tmassWBoson_ < 0)
193  operator()(b, mu, met);
194  return tmassWBoson_;
195 }
double tmassWBoson_
cache of W boson transverse mass estimate
void operator()(const std::vector< pat::Jet > &jets)
double b
Definition: hdecay.h:120

Member Data Documentation

◆ failed_

bool Calculate_miniAOD::failed_
private

indicate failed associations

Definition at line 111 of file TopDQMHelpers.h.

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

◆ massBTopQuark_

double Calculate_miniAOD::massBTopQuark_
private

cache of b-tagged top quark mass estimate

Definition at line 121 of file TopDQMHelpers.h.

Referenced by massBTopQuark(), and operator2().

◆ massTopQuark_

double Calculate_miniAOD::massTopQuark_
private

cache of top quark mass estimate

Definition at line 119 of file TopDQMHelpers.h.

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

◆ massWBoson_

double Calculate_miniAOD::massWBoson_
private

cache of w boson mass estimate

Definition at line 117 of file TopDQMHelpers.h.

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

◆ maxNJets_

int Calculate_miniAOD::maxNJets_
private

max. number of jets to be considered

Definition at line 113 of file TopDQMHelpers.h.

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

◆ tmassTopQuark_

double Calculate_miniAOD::tmassTopQuark_
private

cache of top quark transverse mass estimate

Definition at line 125 of file TopDQMHelpers.h.

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

◆ tmassWBoson_

double Calculate_miniAOD::tmassWBoson_
private

cache of W boson transverse mass estimate

Definition at line 123 of file TopDQMHelpers.h.

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

◆ wMass_

double Calculate_miniAOD::wMass_
private

paramater of the w boson mass

Definition at line 115 of file TopDQMHelpers.h.

Referenced by operator()().