CMS 3D CMS Logo

TtSemiLepHitFitProducer.cc
Go to the documentation of this file.
5 
10 
16 
17 template <typename LeptonCollection>
19 public:
21 
22 private:
23  // produce
24  void produce(edm::Event&, const edm::EventSetup&) override;
25 
29 
31  int maxNJets_;
33  int maxNComb_;
34 
36  double maxEtaMu_;
38  double maxEtaEle_;
40  double maxEtaJet_;
41 
49  bool useBTag_;
50 
52  double mW_;
53  double mTop_;
54 
57 
59  double jes_;
60  double jesB_;
61 
62  struct FitResult {
63  int Status;
64  double Chi2;
65  double Prob;
66  double MT;
67  double SigMT;
74  std::vector<int> JetCombi;
75  bool operator<(const FitResult& rhs) { return Chi2 < rhs.Chi2; };
76  };
77 
79 
86 
91 
92  std::unique_ptr<PatHitFit> HitFit;
93 };
94 
95 template <typename LeptonCollection>
97  : jetsToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
98  lepsToken_(consumes<LeptonCollection>(cfg.getParameter<edm::InputTag>("leps"))),
99  metsToken_(consumes<std::vector<pat::MET> >(cfg.getParameter<edm::InputTag>("mets"))),
100  maxNJets_(cfg.getParameter<int>("maxNJets")),
101  maxNComb_(cfg.getParameter<int>("maxNComb")),
102  bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
103  minBTagValueBJet_(cfg.getParameter<double>("minBDiscBJets")),
104  maxBTagValueNonBJet_(cfg.getParameter<double>("maxBDiscLightJets")),
105  useBTag_(cfg.getParameter<bool>("useBTagging")),
106  mW_(cfg.getParameter<double>("mW")),
107  mTop_(cfg.getParameter<double>("mTop")),
108  jetCorrectionLevel_(cfg.getParameter<std::string>("jetCorrectionLevel")),
109  jes_(cfg.getParameter<double>("jes")),
110  jesB_(cfg.getParameter<double>("jesB")),
111 
112  // The following five initializers read the config parameters for the
113  // ASCII text files which contains the physics object resolutions.
114  hitfitDefault_(cfg.getUntrackedParameter<edm::FileInPath>(
115  std::string("hitfitDefault"),
116  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/setting/RunHitFitConfiguration.txt")))),
117  hitfitElectronResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
118  std::string("hitfitElectronResolution"),
119  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafElectronResolution.txt")))),
120  hitfitMuonResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
121  std::string("hitfitMuonResolution"),
122  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafMuonResolution.txt")))),
123  hitfitUdscJetResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
124  std::string("hitfitUdscJetResolution"),
125  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafUdscJetResolution.txt")))),
126  hitfitBJetResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
127  std::string("hitfitBJetResolution"),
128  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafBJetResolution.txt")))),
129  hitfitMETResolution_(cfg.getUntrackedParameter<edm::FileInPath>(
130  std::string("hitfitMETResolution"),
131  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafKtResolution.txt")))),
132 
133  // The following four initializers instantiate the translator between PAT objects
134  // and HitFit objects using the ASCII text files which contains the resolutions.
135  electronTranslator_(hitfitElectronResolution_.fullPath()),
136  muonTranslator_(hitfitMuonResolution_.fullPath()),
137  jetTranslator_(
138  hitfitUdscJetResolution_.fullPath(), hitfitBJetResolution_.fullPath(), jetCorrectionLevel_, jes_, jesB_),
139  metTranslator_(hitfitMETResolution_.fullPath())
140 
141 {
142  // Create an instance of RunHitFit and initialize it.
143  HitFit = std::make_unique<PatHitFit>(
145 
146  maxEtaMu_ = 2.4;
147  maxEtaEle_ = 2.5;
148  maxEtaJet_ = 2.5;
149 
150  edm::LogVerbatim("TopHitFit") << "\n"
151  << "+++++++++++ TtSemiLepHitFitProducer ++++++++++++ \n"
152  << " Due to the eta ranges for which resolutions \n"
153  << " are provided in \n"
154  << " TopQuarkAnalysis/TopHitFit/data/resolution/ \n"
155  << " so far, the following cuts are currently \n"
156  << " implemented in the TtSemiLepHitFitProducer: \n"
157  << " |eta(muons )| <= " << maxEtaMu_ << " \n"
158  << " |eta(electrons)| <= " << maxEtaEle_ << " \n"
159  << " |eta(jets )| <= " << maxEtaJet_ << " \n"
160  << "++++++++++++++++++++++++++++++++++++++++++++++++ \n";
161 
162  produces<std::vector<pat::Particle> >("PartonsHadP");
163  produces<std::vector<pat::Particle> >("PartonsHadQ");
164  produces<std::vector<pat::Particle> >("PartonsHadB");
165  produces<std::vector<pat::Particle> >("PartonsLepB");
166  produces<std::vector<pat::Particle> >("Leptons");
167  produces<std::vector<pat::Particle> >("Neutrinos");
168 
169  produces<std::vector<std::vector<int> > >();
170  produces<std::vector<double> >("Chi2");
171  produces<std::vector<double> >("Prob");
172  produces<std::vector<double> >("MT");
173  produces<std::vector<double> >("SigMT");
174  produces<std::vector<int> >("Status");
175  produces<int>("NumberOfConsideredJets");
176 }
177 
178 template <typename LeptonCollection>
180  std::unique_ptr<std::vector<pat::Particle> > pPartonsHadP(new std::vector<pat::Particle>);
181  std::unique_ptr<std::vector<pat::Particle> > pPartonsHadQ(new std::vector<pat::Particle>);
182  std::unique_ptr<std::vector<pat::Particle> > pPartonsHadB(new std::vector<pat::Particle>);
183  std::unique_ptr<std::vector<pat::Particle> > pPartonsLepB(new std::vector<pat::Particle>);
184  std::unique_ptr<std::vector<pat::Particle> > pLeptons(new std::vector<pat::Particle>);
185  std::unique_ptr<std::vector<pat::Particle> > pNeutrinos(new std::vector<pat::Particle>);
186 
187  std::unique_ptr<std::vector<std::vector<int> > > pCombi(new std::vector<std::vector<int> >);
188  std::unique_ptr<std::vector<double> > pChi2(new std::vector<double>);
189  std::unique_ptr<std::vector<double> > pProb(new std::vector<double>);
190  std::unique_ptr<std::vector<double> > pMT(new std::vector<double>);
191  std::unique_ptr<std::vector<double> > pSigMT(new std::vector<double>);
192  std::unique_ptr<std::vector<int> > pStatus(new std::vector<int>);
193  std::unique_ptr<int> pJetsConsidered(new int);
194 
196  evt.getByToken(jetsToken_, jets);
197 
199  evt.getByToken(metsToken_, mets);
200 
202  evt.getByToken(lepsToken_, leps);
203 
204  // -----------------------------------------------------
205  // skip events with no appropriate lepton candidate in
206  // or empty MET or less jets than partons
207  // -----------------------------------------------------
208 
209  const unsigned int nPartons = 4;
210 
211  // Clear the internal state
212  HitFit->clear();
213 
214  // Add lepton into HitFit
215  bool foundLepton = false;
216  if (!leps->empty()) {
217  double maxEtaLep = maxEtaMu_;
218  if (!dynamic_cast<const reco::Muon*>(&((*leps)[0]))) // assume electron if it is not a muon
219  maxEtaLep = maxEtaEle_;
220  for (unsigned iLep = 0; iLep < (*leps).size() && !foundLepton; ++iLep) {
221  if (std::abs((*leps)[iLep].eta()) <= maxEtaLep) {
222  HitFit->AddLepton((*leps)[iLep]);
223  foundLepton = true;
224  }
225  }
226  }
227 
228  // Add jets into HitFit
229  unsigned int nJetsFound = 0;
230  for (unsigned iJet = 0; iJet < (*jets).size() && (int)nJetsFound != maxNJets_; ++iJet) {
231  double jet_eta = (*jets)[iJet].eta();
232  if ((*jets)[iJet].isCaloJet()) {
233  jet_eta = ((reco::CaloJet*)((*jets)[iJet]).originalObject())->detectorP4().eta();
234  }
235  if (std::abs(jet_eta) <= maxEtaJet_) {
236  HitFit->AddJet((*jets)[iJet]);
237  nJetsFound++;
238  }
239  }
240  *pJetsConsidered = nJetsFound;
241 
242  // Add missing transverse energy into HitFit
243  if (!mets->empty())
244  HitFit->SetMet((*mets)[0]);
245 
246  if (!foundLepton || mets->empty() || nJetsFound < nPartons) {
247  // the kinFit getters return empty objects here
248  pPartonsHadP->push_back(pat::Particle());
249  pPartonsHadQ->push_back(pat::Particle());
250  pPartonsHadB->push_back(pat::Particle());
251  pPartonsLepB->push_back(pat::Particle());
252  pLeptons->push_back(pat::Particle());
253  pNeutrinos->push_back(pat::Particle());
254  // indices referring to the jet combination
255  std::vector<int> invalidCombi;
256  for (unsigned int i = 0; i < nPartons; ++i)
257  invalidCombi.push_back(-1);
258  pCombi->push_back(invalidCombi);
259  // chi2
260  pChi2->push_back(-1.);
261  // chi2 probability
262  pProb->push_back(-1.);
263  // fitted top mass
264  pMT->push_back(-1.);
265  pSigMT->push_back(-1.);
266  // status of the fitter
267  pStatus->push_back(-1);
268  // feed out all products
269  evt.put(std::move(pCombi));
270  evt.put(std::move(pPartonsHadP), "PartonsHadP");
271  evt.put(std::move(pPartonsHadQ), "PartonsHadQ");
272  evt.put(std::move(pPartonsHadB), "PartonsHadB");
273  evt.put(std::move(pPartonsLepB), "PartonsLepB");
274  evt.put(std::move(pLeptons), "Leptons");
275  evt.put(std::move(pNeutrinos), "Neutrinos");
276  evt.put(std::move(pChi2), "Chi2");
277  evt.put(std::move(pProb), "Prob");
278  evt.put(std::move(pMT), "MT");
279  evt.put(std::move(pSigMT), "SigMT");
280  evt.put(std::move(pStatus), "Status");
281  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
282  return;
283  }
284 
285  std::list<FitResult> FitResultList;
286 
287  //
288  // BEGIN DECLARATION OF VARIABLES FROM KINEMATIC FIT
289  //
290 
291  // In this part are variables from the
292  // kinematic fit procedure
293 
294  // Number of all permutations of the event
295  size_t nHitFit = 0;
296 
297  // Number of jets in the event
298  size_t nHitFitJet = 0;
299 
300  // Results of the fit for all jet permutations of the event
301  std::vector<hitfit::Fit_Result> hitFitResult;
302 
303  //
304  // R U N H I T F I T
305  //
306  // Run the kinematic fit and get how many permutations are possible
307  // in the fit
308 
309  nHitFit = HitFit->FitAllPermutation();
310 
311  //
312  // BEGIN PART WHICH EXTRACTS INFORMATION FROM HITFIT
313  //
314 
315  // Get the number of jets
316  nHitFitJet = HitFit->GetUnfittedEvent()[0].njets();
317 
318  // Get the fit results for all permutations
319  hitFitResult = HitFit->GetFitAllPermutation();
320 
321  // Loop over all permutations and extract the information
322  for (size_t fit = 0; fit != nHitFit; ++fit) {
323  // Get the event after the fit
324  hitfit::Lepjets_Event fittedEvent = hitFitResult[fit].ev();
325 
326  /*
327  Get jet permutation according to TQAF convention
328  11 : leptonic b
329  12 : hadronic b
330  13 : hadronic W
331  14 : hadronic W
332  */
333  std::vector<int> hitCombi(4);
334  for (size_t jet = 0; jet != nHitFitJet; ++jet) {
335  int jet_type = fittedEvent.jet(jet).type();
336 
337  switch (jet_type) {
338  case 11:
339  hitCombi[TtSemiLepEvtPartons::LepB] = jet;
340  break;
341  case 12:
342  hitCombi[TtSemiLepEvtPartons::HadB] = jet;
343  break;
344  case 13:
345  hitCombi[TtSemiLepEvtPartons::LightQ] = jet;
346  break;
347  case 14:
349  break;
350  }
351  }
352 
353  // Store the kinematic quantities in the corresponding containers.
354 
355  hitfit::Lepjets_Event_Jet hadP_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LightQ]);
356  hitfit::Lepjets_Event_Jet hadQ_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LightQBar]);
357  hitfit::Lepjets_Event_Jet hadB_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::HadB]);
358  hitfit::Lepjets_Event_Jet lepB_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LepB]);
359  hitfit::Lepjets_Event_Lep lepL_ = fittedEvent.lep(0);
360 
361  /*
363  std::string bTagAlgo_;
365  double minBTagValueBJet_;
367  double maxBTagValueNonBJet_;
369  bool useBTag_;
370  */
371 
372  if (hitFitResult[fit].chisq() > 0 // only take into account converged fits
373  && (!useBTag_ ||
374  (useBTag_ // use btag information if chosen
375  && jets->at(hitCombi[TtSemiLepEvtPartons::LightQ]).bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_ &&
376  jets->at(hitCombi[TtSemiLepEvtPartons::LightQBar]).bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_ &&
377  jets->at(hitCombi[TtSemiLepEvtPartons::HadB]).bDiscriminator(bTagAlgo_) > minBTagValueBJet_ &&
378  jets->at(hitCombi[TtSemiLepEvtPartons::LepB]).bDiscriminator(bTagAlgo_) > minBTagValueBJet_))) {
380  result.Status = 0;
381  result.Chi2 = hitFitResult[fit].chisq();
382  result.Prob = exp(-1.0 * (hitFitResult[fit].chisq()) / 2.0);
383  result.MT = hitFitResult[fit].mt();
384  result.SigMT = hitFitResult[fit].sigmt();
386  0, math::XYZTLorentzVector(hadB_.p().x(), hadB_.p().y(), hadB_.p().z(), hadB_.p().t()), math::XYZPoint()));
388  0, math::XYZTLorentzVector(hadP_.p().x(), hadP_.p().y(), hadP_.p().z(), hadP_.p().t()), math::XYZPoint()));
390  0, math::XYZTLorentzVector(hadQ_.p().x(), hadQ_.p().y(), hadQ_.p().z(), hadQ_.p().t()), math::XYZPoint()));
392  0, math::XYZTLorentzVector(lepB_.p().x(), lepB_.p().y(), lepB_.p().z(), lepB_.p().t()), math::XYZPoint()));
394  0, math::XYZTLorentzVector(lepL_.p().x(), lepL_.p().y(), lepL_.p().z(), lepL_.p().t()), math::XYZPoint()));
396  0,
398  fittedEvent.met().x(), fittedEvent.met().y(), fittedEvent.met().z(), fittedEvent.met().t()),
399  math::XYZPoint()));
400  result.JetCombi = hitCombi;
401 
402  FitResultList.push_back(result);
403  }
404  }
405 
406  // sort results w.r.t. chi2 values
407  FitResultList.sort();
408 
409  // -----------------------------------------------------
410  // feed out result
411  // starting with the JetComb having the smallest chi2
412  // -----------------------------------------------------
413 
414  if (((unsigned)FitResultList.size()) < 1) { // in case no fit results were stored in the list (all fits aborted)
415  pPartonsHadP->push_back(pat::Particle());
416  pPartonsHadQ->push_back(pat::Particle());
417  pPartonsHadB->push_back(pat::Particle());
418  pPartonsLepB->push_back(pat::Particle());
419  pLeptons->push_back(pat::Particle());
420  pNeutrinos->push_back(pat::Particle());
421  // indices referring to the jet combination
422  std::vector<int> invalidCombi;
423  for (unsigned int i = 0; i < nPartons; ++i)
424  invalidCombi.push_back(-1);
425  pCombi->push_back(invalidCombi);
426  // chi2
427  pChi2->push_back(-1.);
428  // chi2 probability
429  pProb->push_back(-1.);
430  // fitted top mass
431  pMT->push_back(-1.);
432  pSigMT->push_back(-1.);
433  // status of the fitter
434  pStatus->push_back(-1);
435  } else {
436  unsigned int iComb = 0;
437  for (typename std::list<FitResult>::const_iterator result = FitResultList.begin(); result != FitResultList.end();
438  ++result) {
439  if (maxNComb_ >= 1 && iComb == (unsigned int)maxNComb_)
440  break;
441  iComb++;
442  // partons
443  pPartonsHadP->push_back(result->HadP);
444  pPartonsHadQ->push_back(result->HadQ);
445  pPartonsHadB->push_back(result->HadB);
446  pPartonsLepB->push_back(result->LepB);
447  // lepton
448  pLeptons->push_back(result->LepL);
449  // neutrino
450  pNeutrinos->push_back(result->LepN);
451  // indices referring to the jet combination
452  pCombi->push_back(result->JetCombi);
453  // chi2
454  pChi2->push_back(result->Chi2);
455  // chi2 probability
456  pProb->push_back(result->Prob);
457  // fitted top mass
458  pMT->push_back(result->MT);
459  pSigMT->push_back(result->SigMT);
460  // status of the fitter
461  pStatus->push_back(result->Status);
462  }
463  }
464  evt.put(std::move(pCombi));
465  evt.put(std::move(pPartonsHadP), "PartonsHadP");
466  evt.put(std::move(pPartonsHadQ), "PartonsHadQ");
467  evt.put(std::move(pPartonsHadB), "PartonsHadB");
468  evt.put(std::move(pPartonsLepB), "PartonsLepB");
469  evt.put(std::move(pLeptons), "Leptons");
470  evt.put(std::move(pNeutrinos), "Neutrinos");
471  evt.put(std::move(pChi2), "Chi2");
472  evt.put(std::move(pProb), "Prob");
473  evt.put(std::move(pMT), "MT");
474  evt.put(std::move(pSigMT), "SigMT");
475  evt.put(std::move(pStatus), "Status");
476  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
477 }
478 
483 
Log< level::Info, true > LogVerbatim
int maxNComb_
maximal number of combinations to be written to the event
double maxEtaEle_
maximum eta value for electrons, needed to limited range in which resolutions are provided ...
double maxBTagValueNonBJet_
max value of bTag for a non-b-jet
static const unsigned int nPartons
hitfit::JetTranslatorBase< pat::Jet > jetTranslator_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
Jets made from CaloTowers.
Definition: CaloJet.h:27
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information: ...
std::string fullPath() const
Definition: FileInPath.cc:161
hitfit::LeptonTranslatorBase< pat::Electron > electronTranslator_
TtSemiLepHitFitProducer(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
hitfit::LeptonTranslatorBase< pat::Muon > muonTranslator_
Lepjets_Event_Jet & jet(std::vector< Lepjets_Event_Jet >::size_type i)
Return a reference to jet at index position i.
double maxEtaJet_
maximum eta value for jets, needed to limited range in which resolutions are provided ...
Definition: HeavyIon.h:7
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int & type()
Return a reference to the type code.
edm::FileInPath hitfitUdscJetResolution_
void produce(edm::Event &, const edm::EventSetup &) override
Definition: Jet.py:1
A class to hold functions to calculate kinematic quantities of interest in events.
bool useBTag_
switch to tell whether to use b-tagging or not
TtSemiLepHitFitProducer< std::vector< pat::Electron > > TtSemiLepHitFitProducerElectron
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:62
std::unique_ptr< PatHitFit > HitFit
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Template class of function object to translate missing transverse energy object to HitFit&#39;s Fourvec o...
A class to represent a jet in an instance of Lepjets_Event class. The class is derived from the Lepje...
int maxNJets_
maximal number of jets (-1 possible to indicate &#39;all&#39;)
Template class of experiment-independent interface to HitFit.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
Template class of experiment-independent interface to HitFit. This class is intended to be used insid...
Definition: RunHitFit.h:153
hitfit::METTranslatorBase< pat::MET > metTranslator_
Template class of function object to translate jet physics object to HitFit&#39;s Lepjets_Event_Jet objec...
Fourvec & met()
Return a reference to the missing transverse energy.
Lepjets_Event_Lep & lep(std::vector< Lepjets_Event_Lep >::size_type i)
Return a reference to lepton at index position i.
edm::EDGetTokenT< LeptonCollection > lepsToken_
Template class of function object to translate lepton physics object to HitFit&#39;s Lepjets_Event_Lep ob...
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::FileInPath hitfitElectronResolution_
Analysis-level particle class.
Definition: Particle.h:30
std::string jetCorrectionLevel_
jet correction level
Fourvec & p()
Return a reference to the four-momentum.
TtSemiLepHitFitProducer< std::vector< pat::Muon > > TtSemiLepHitFitProducerMuon
hitfit::RunHitFit< pat::Electron, pat::Muon, pat::Jet, pat::MET > PatHitFit
HLT enums.
double maxEtaMu_
maximum eta value for muons, needed to limited range in which resolutions are provided ...
Definition: Chi2.h:15
double minBTagValueBJet_
min value of bTag for a b-jet
std::string bTagAlgo_
input tag for b-tagging algorithm
def move(src, dest)
Definition: eostools.py:511
double jes_
jet energy scale