CMS 3D CMS Logo

TtSemiLepKinFitProducer.cc
Go to the documentation of this file.
4 
8 
9 template <typename LeptonCollection>
11 public:
13 
14 private:
15  // produce
16  void produce(edm::Event&, const edm::EventSetup&) override;
17 
18  // convert unsigned to Param
20  // convert unsigned to Param
22  // convert unsigned to Param
23  std::vector<TtSemiLepKinFitter::Constraint> constraints(std::vector<unsigned>&);
24  // helper function for b-tagging
25  bool doBTagging(bool& useBTag_,
26  const edm::Handle<std::vector<pat::Jet>>& jets,
27  const std::vector<int>& combi,
29  double& minBTagValueBJets_,
30  double& maxBTagValueNonBJets_);
31 
35 
46  bool useBTag_;
48  int maxNJets_;
50  int maxNComb_;
51 
53  unsigned int maxNrIter_;
55  double maxDeltaS_;
57  double maxF_;
58  unsigned int jetParam_;
59  unsigned int lepParam_;
60  unsigned int metParam_;
62  std::vector<unsigned> constraints_;
63  double mW_;
64  double mTop_;
66  std::vector<double> jetEnergyResolutionScaleFactors_;
67  std::vector<double> jetEnergyResolutionEtaBinning_;
69  std::vector<edm::ParameterSet> udscResolutions_;
70  std::vector<edm::ParameterSet> bResolutions_;
71  std::vector<edm::ParameterSet> lepResolutions_;
72  std::vector<edm::ParameterSet> metResolutions_;
73 
74  std::unique_ptr<TtSemiLepKinFitter> fitter;
75 
76  struct KinFitResult {
77  int Status;
78  double Chi2;
79  double Prob;
86  std::vector<int> JetCombi;
87  bool operator<(const KinFitResult& rhs) { return Chi2 < rhs.Chi2; };
88  };
89 };
90 
91 template <typename LeptonCollection>
93  : jetsToken_(consumes<std::vector<pat::Jet>>(cfg.getParameter<edm::InputTag>("jets"))),
94  lepsToken_(consumes<LeptonCollection>(cfg.getParameter<edm::InputTag>("leps"))),
95  metsToken_(consumes<std::vector<pat::MET>>(cfg.getParameter<edm::InputTag>("mets"))),
96  matchToken_(mayConsume<std::vector<std::vector<int>>>(cfg.getParameter<edm::InputTag>("match"))),
97  useOnlyMatch_(cfg.getParameter<bool>("useOnlyMatch")),
98  bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
99  minBTagValueBJet_(cfg.getParameter<double>("minBDiscBJets")),
100  maxBTagValueNonBJet_(cfg.getParameter<double>("maxBDiscLightJets")),
101  useBTag_(cfg.getParameter<bool>("useBTagging")),
102  maxNJets_(cfg.getParameter<int>("maxNJets")),
103  maxNComb_(cfg.getParameter<int>("maxNComb")),
104  maxNrIter_(cfg.getParameter<unsigned>("maxNrIter")),
105  maxDeltaS_(cfg.getParameter<double>("maxDeltaS")),
106  maxF_(cfg.getParameter<double>("maxF")),
107  jetParam_(cfg.getParameter<unsigned>("jetParametrisation")),
108  lepParam_(cfg.getParameter<unsigned>("lepParametrisation")),
109  metParam_(cfg.getParameter<unsigned>("metParametrisation")),
110  constraints_(cfg.getParameter<std::vector<unsigned>>("constraints")),
111  mW_(cfg.getParameter<double>("mW")),
112  mTop_(cfg.getParameter<double>("mTop")),
113  jetEnergyResolutionScaleFactors_(cfg.getParameter<std::vector<double>>("jetEnergyResolutionScaleFactors")),
114  jetEnergyResolutionEtaBinning_(cfg.getParameter<std::vector<double>>("jetEnergyResolutionEtaBinning")),
115  udscResolutions_(0),
116  bResolutions_(0),
117  lepResolutions_(0),
118  metResolutions_(0) {
119  if (cfg.exists("udscResolutions") && cfg.exists("bResolutions") && cfg.exists("lepResolutions") &&
120  cfg.exists("metResolutions")) {
121  udscResolutions_ = cfg.getParameter<std::vector<edm::ParameterSet>>("udscResolutions");
122  bResolutions_ = cfg.getParameter<std::vector<edm::ParameterSet>>("bResolutions");
123  lepResolutions_ = cfg.getParameter<std::vector<edm::ParameterSet>>("lepResolutions");
124  metResolutions_ = cfg.getParameter<std::vector<edm::ParameterSet>>("metResolutions");
125  } else if (cfg.exists("udscResolutions") || cfg.exists("bResolutions") || cfg.exists("lepResolutions") ||
126  cfg.exists("metResolutions")) {
127  throw cms::Exception("Configuration") << "Parameters 'udscResolutions', 'bResolutions', 'lepResolutions', "
128  "'metResolutions' should be used together.\n";
129  }
130 
131  fitter = std::make_unique<TtSemiLepKinFitter>(param(jetParam_),
132  param(lepParam_),
133  param(metParam_),
134  maxNrIter_,
135  maxDeltaS_,
136  maxF_,
138  mW_,
139  mTop_,
141  &bResolutions_,
146 
147  produces<std::vector<pat::Particle>>("PartonsHadP");
148  produces<std::vector<pat::Particle>>("PartonsHadQ");
149  produces<std::vector<pat::Particle>>("PartonsHadB");
150  produces<std::vector<pat::Particle>>("PartonsLepB");
151  produces<std::vector<pat::Particle>>("Leptons");
152  produces<std::vector<pat::Particle>>("Neutrinos");
153 
154  produces<std::vector<std::vector<int>>>();
155  produces<std::vector<double>>("Chi2");
156  produces<std::vector<double>>("Prob");
157  produces<std::vector<int>>("Status");
158 
159  produces<int>("NumberOfConsideredJets");
160 }
161 
162 template <typename LeptonCollection>
164  const edm::Handle<std::vector<pat::Jet>>& jets,
165  const std::vector<int>& combi,
166  std::string& bTagAlgo_,
167  double& minBTagValueBJet_,
168  double& maxBTagValueNonBJet_) {
169  if (!useBTag_) {
170  return true;
171  }
172  if (useBTag_ && (*jets)[combi[TtSemiLepEvtPartons::HadB]].bDiscriminator(bTagAlgo_) >= minBTagValueBJet_ &&
173  (*jets)[combi[TtSemiLepEvtPartons::LepB]].bDiscriminator(bTagAlgo_) >= minBTagValueBJet_ &&
174  (*jets)[combi[TtSemiLepEvtPartons::LightQ]].bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_ &&
175  (*jets)[combi[TtSemiLepEvtPartons::LightQBar]].bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_) {
176  return true;
177  } else {
178  return false;
179  }
180 }
181 
182 template <typename LeptonCollection>
184  std::unique_ptr<std::vector<pat::Particle>> pPartonsHadP(new std::vector<pat::Particle>);
185  std::unique_ptr<std::vector<pat::Particle>> pPartonsHadQ(new std::vector<pat::Particle>);
186  std::unique_ptr<std::vector<pat::Particle>> pPartonsHadB(new std::vector<pat::Particle>);
187  std::unique_ptr<std::vector<pat::Particle>> pPartonsLepB(new std::vector<pat::Particle>);
188  std::unique_ptr<std::vector<pat::Particle>> pLeptons(new std::vector<pat::Particle>);
189  std::unique_ptr<std::vector<pat::Particle>> pNeutrinos(new std::vector<pat::Particle>);
190 
191  std::unique_ptr<std::vector<std::vector<int>>> pCombi(new std::vector<std::vector<int>>);
192  std::unique_ptr<std::vector<double>> pChi2(new std::vector<double>);
193  std::unique_ptr<std::vector<double>> pProb(new std::vector<double>);
194  std::unique_ptr<std::vector<int>> pStatus(new std::vector<int>);
195 
196  std::unique_ptr<int> pJetsConsidered(new int);
197 
198  const edm::Handle<std::vector<pat::Jet>>& jets = evt.getHandle(jetsToken_);
199 
200  const edm::Handle<std::vector<pat::MET>>& mets = evt.getHandle(metsToken_);
201 
202  const edm::Handle<LeptonCollection>& leps = evt.getHandle(lepsToken_);
203 
204  const unsigned int nPartons = 4;
205 
206  std::vector<int> match;
207  bool invalidMatch = false;
208  if (useOnlyMatch_) {
209  *pJetsConsidered = nPartons;
210  const edm::Handle<std::vector<std::vector<int>>>& matchHandle = evt.getHandle(matchToken_);
211  match = *(matchHandle->begin());
212  // check if match is valid
213  if (match.size() != nPartons)
214  invalidMatch = true;
215  else {
216  for (unsigned int idx = 0; idx < match.size(); ++idx) {
217  if (match[idx] < 0 || match[idx] >= (int)jets->size()) {
218  invalidMatch = true;
219  break;
220  }
221  }
222  }
223  }
224 
225  // -----------------------------------------------------
226  // skip events with no appropriate lepton candidate in
227  // or empty MET or less jets than partons or invalid match
228  // -----------------------------------------------------
229 
230  if (leps->empty() || mets->empty() || jets->size() < nPartons || invalidMatch) {
231  // the kinFit getters return empty objects here
232  pPartonsHadP->push_back(fitter->fittedHadP());
233  pPartonsHadQ->push_back(fitter->fittedHadQ());
234  pPartonsHadB->push_back(fitter->fittedHadB());
235  pPartonsLepB->push_back(fitter->fittedLepB());
236  pLeptons->push_back(fitter->fittedLepton());
237  pNeutrinos->push_back(fitter->fittedNeutrino());
238  // indices referring to the jet combination
239  std::vector<int> invalidCombi;
240  for (unsigned int i = 0; i < nPartons; ++i)
241  invalidCombi.push_back(-1);
242  pCombi->push_back(invalidCombi);
243  // chi2
244  pChi2->push_back(-1.);
245  // chi2 probability
246  pProb->push_back(-1.);
247  // status of the fitter
248  pStatus->push_back(-1);
249  // number of jets
250  *pJetsConsidered = jets->size();
251  // feed out all products
252  evt.put(std::move(pCombi));
253  evt.put(std::move(pPartonsHadP), "PartonsHadP");
254  evt.put(std::move(pPartonsHadQ), "PartonsHadQ");
255  evt.put(std::move(pPartonsHadB), "PartonsHadB");
256  evt.put(std::move(pPartonsLepB), "PartonsLepB");
257  evt.put(std::move(pLeptons), "Leptons");
258  evt.put(std::move(pNeutrinos), "Neutrinos");
259  evt.put(std::move(pChi2), "Chi2");
260  evt.put(std::move(pProb), "Prob");
261  evt.put(std::move(pStatus), "Status");
262  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
263  return;
264  }
265 
266  // -----------------------------------------------------
267  // analyze different jet combinations using the KinFitter
268  // (or only a given jet combination if useOnlyMatch=true)
269  // -----------------------------------------------------
270 
271  std::vector<int> jetIndices;
272  if (!useOnlyMatch_) {
273  for (unsigned int i = 0; i < jets->size(); ++i) {
274  if (maxNJets_ >= (int)nPartons && maxNJets_ == (int)i) {
275  *pJetsConsidered = i;
276  break;
277  }
278  jetIndices.push_back(i);
279  }
280  }
281 
282  std::vector<int> combi;
283  for (unsigned int i = 0; i < nPartons; ++i) {
284  if (useOnlyMatch_)
285  combi.push_back(match[i]);
286  else
287  combi.push_back(i);
288  }
289 
290  std::list<KinFitResult> FitResultList;
291 
292  do {
293  for (int cnt = 0; cnt < TMath::Factorial(combi.size()); ++cnt) {
294  // take into account indistinguishability of the two jets from the hadr. W decay,
295  // reduces combinatorics by a factor of 2
296  if ((combi[TtSemiLepEvtPartons::LightQ] < combi[TtSemiLepEvtPartons::LightQBar] || useOnlyMatch_) &&
297  doBTagging(useBTag_, jets, combi, bTagAlgo_, minBTagValueBJet_, maxBTagValueNonBJet_)) {
298  std::vector<pat::Jet> jetCombi;
299  jetCombi.resize(nPartons);
300  jetCombi[TtSemiLepEvtPartons::LightQ] = (*jets)[combi[TtSemiLepEvtPartons::LightQ]];
302  jetCombi[TtSemiLepEvtPartons::HadB] = (*jets)[combi[TtSemiLepEvtPartons::HadB]];
303  jetCombi[TtSemiLepEvtPartons::LepB] = (*jets)[combi[TtSemiLepEvtPartons::LepB]];
304 
305  // do the kinematic fit
306  const int status = fitter->fit(jetCombi, (*leps)[0], (*mets)[0]);
307 
308  if (status == 0) { // only take into account converged fits
310  result.Status = status;
311  result.Chi2 = fitter->fitS();
312  result.Prob = fitter->fitProb();
313  result.HadB = fitter->fittedHadB();
314  result.HadP = fitter->fittedHadP();
315  result.HadQ = fitter->fittedHadQ();
316  result.LepB = fitter->fittedLepB();
317  result.LepL = fitter->fittedLepton();
318  result.LepN = fitter->fittedNeutrino();
319  result.JetCombi = combi;
320 
321  FitResultList.push_back(result);
322  }
323  }
324  if (useOnlyMatch_)
325  break; // don't go through combinatorics if useOnlyMatch was chosen
326  next_permutation(combi.begin(), combi.end());
327  }
328  if (useOnlyMatch_)
329  break; // don't go through combinatorics if useOnlyMatch was chosen
330  } while (stdcomb::next_combination(jetIndices.begin(), jetIndices.end(), combi.begin(), combi.end()));
331 
332  // sort results w.r.t. chi2 values
333  FitResultList.sort();
334 
335  // -----------------------------------------------------
336  // feed out result
337  // starting with the JetComb having the smallest chi2
338  // -----------------------------------------------------
339 
340  if ((unsigned)FitResultList.size() < 1) { // in case no fit results were stored in the list (all fits aborted)
341  pPartonsHadP->push_back(fitter->fittedHadP());
342  pPartonsHadQ->push_back(fitter->fittedHadQ());
343  pPartonsHadB->push_back(fitter->fittedHadB());
344  pPartonsLepB->push_back(fitter->fittedLepB());
345  pLeptons->push_back(fitter->fittedLepton());
346  pNeutrinos->push_back(fitter->fittedNeutrino());
347  // indices referring to the jet combination
348  std::vector<int> invalidCombi;
349  for (unsigned int i = 0; i < nPartons; ++i)
350  invalidCombi.push_back(-1);
351  pCombi->push_back(invalidCombi);
352  // chi2
353  pChi2->push_back(-1.);
354  // chi2 probability
355  pProb->push_back(-1.);
356  // status of the fitter
357  pStatus->push_back(-1);
358  } else {
359  unsigned int iComb = 0;
360  for (typename std::list<KinFitResult>::const_iterator result = FitResultList.begin(); result != FitResultList.end();
361  ++result) {
362  if (maxNComb_ >= 1 && iComb == (unsigned int)maxNComb_)
363  break;
364  iComb++;
365  // partons
366  pPartonsHadP->push_back(result->HadP);
367  pPartonsHadQ->push_back(result->HadQ);
368  pPartonsHadB->push_back(result->HadB);
369  pPartonsLepB->push_back(result->LepB);
370  // lepton
371  pLeptons->push_back(result->LepL);
372  // neutrino
373  pNeutrinos->push_back(result->LepN);
374  // indices referring to the jet combination
375  pCombi->push_back(result->JetCombi);
376  // chi2
377  pChi2->push_back(result->Chi2);
378  // chi2 probability
379  pProb->push_back(result->Prob);
380  // status of the fitter
381  pStatus->push_back(result->Status);
382  }
383  }
384  evt.put(std::move(pCombi));
385  evt.put(std::move(pPartonsHadP), "PartonsHadP");
386  evt.put(std::move(pPartonsHadQ), "PartonsHadQ");
387  evt.put(std::move(pPartonsHadB), "PartonsHadB");
388  evt.put(std::move(pPartonsLepB), "PartonsLepB");
389  evt.put(std::move(pLeptons), "Leptons");
390  evt.put(std::move(pNeutrinos), "Neutrinos");
391  evt.put(std::move(pChi2), "Chi2");
392  evt.put(std::move(pProb), "Prob");
393  evt.put(std::move(pStatus), "Status");
394  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
395 }
396 
397 template <typename LeptonCollection>
400  switch (val) {
403  break;
406  break;
409  break;
410  default:
411  throw cms::Exception("Configuration") << "Chosen jet parametrization is not supported: " << val << "\n";
412  break;
413  }
414  return result;
415 }
416 
417 template <typename LeptonCollection>
420  switch (val) {
423  break;
426  break;
429  break;
432  break;
435  break;
438  break;
441  break;
442  default:
443  throw cms::Exception("Configuration") << "Chosen fit constraint is not supported: " << val << "\n";
444  break;
445  }
446  return result;
447 }
448 
449 template <typename LeptonCollection>
450 std::vector<TtSemiLepKinFitter::Constraint> TtSemiLepKinFitProducer<LeptonCollection>::constraints(
451  std::vector<unsigned>& val) {
452  std::vector<TtSemiLepKinFitter::Constraint> result;
453  for (unsigned i = 0; i < val.size(); ++i) {
454  result.push_back(constraint(val[i]));
455  }
456  return result;
457 }
458 
463 
std::vector< edm::ParameterSet > bResolutions_
std::vector< double > jetEnergyResolutionScaleFactors_
scale factors for jet energy resolution
double maxBTagValueNonBJet_
max value of bTag for a non-b-jet
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
std::string bTagAlgo_
input tag for b-tagging algorithm
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
static const unsigned int nPartons
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
double minBTagValueBJet_
min value of bTag for a b-jet
Param
supported parameterizations
Definition: TopKinFitter.h:22
std::vector< TtSemiLepKinFitter::Constraint > constraints(std::vector< unsigned > &)
bool doBTagging(bool &useBTag_, const edm::Handle< std::vector< pat::Jet >> &jets, const std::vector< int > &combi, std::string &bTagAlgo_, double &minBTagValueBJets_, double &maxBTagValueNonBJets_)
int maxNJets_
maximal number of jets (-1 possible to indicate &#39;all&#39;)
Definition: HeavyIon.h:7
std::unique_ptr< TtSemiLepKinFitter > fitter
std::vector< double > jetEnergyResolutionEtaBinning_
Definition: Jet.py:1
edm::EDGetTokenT< std::vector< std::vector< int > > > matchToken_
std::vector< unsigned > constraints_
constrains
unsigned int maxNrIter_
maximal number of iterations to be performed for the fit
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< edm::ParameterSet > udscResolutions_
config-file-based object resolutions
int maxNComb_
maximal number of combinations to be written to the event
double maxF_
maximal deviation for contstraints
edm::EDGetTokenT< LeptonCollection > lepsToken_
std::vector< edm::ParameterSet > metResolutions_
Analysis-level particle class.
Definition: Particle.h:30
void produce(edm::Event &, const edm::EventSetup &) override
TtSemiLepKinFitter::Constraint constraint(unsigned)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:550
HLT enums.
bool useOnlyMatch_
switch to use only a combination given by another hypothesis
bool useBTag_
switch to tell whether to use b-tagging or not
TtSemiLepKinFitter::Param param(unsigned)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
TtSemiLepKinFitProducer(const edm::ParameterSet &)
double maxDeltaS_
maximal chi2 equivalent
Definition: Chi2.h:15
bool next_combination(BidIt n_begin, BidIt n_end, BidIt r_begin, BidIt r_end)
Definition: combination.h:19
std::vector< edm::ParameterSet > lepResolutions_
def move(src, dest)
Definition: eostools.py:511
Constraint
supported constraints