CMS 3D CMS Logo

TtSemiLepJetCombMaxSumPtWMass.cc
Go to the documentation of this file.
4 
7 
9 public:
11 
12 private:
13  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& setup) const override;
14 
15  bool isValid(const int& idx, const std::vector<pat::Jet>& jets) const {
16  return (0 <= idx && idx < (int)jets.size());
17  };
18 
21  int maxNJets_;
22  double wMass_;
27 };
28 
30  : jetsToken_(consumes<std::vector<pat::Jet>>(cfg.getParameter<edm::InputTag>("jets"))),
31  lepsToken_(consumes<edm::View<reco::RecoCandidate>>(cfg.getParameter<edm::InputTag>("leps"))),
32  maxNJets_(cfg.getParameter<int>("maxNJets")),
33  wMass_(cfg.getParameter<double>("wMass")),
34  useBTagging_(cfg.getParameter<bool>("useBTagging")),
35  bTagAlgorithm_(cfg.getParameter<std::string>("bTagAlgorithm")),
36  minBDiscBJets_(cfg.getParameter<double>("minBDiscBJets")),
37  maxBDiscLightJets_(cfg.getParameter<double>("maxBDiscLightJets")) {
38  if (maxNJets_ < 4 && maxNJets_ != -1)
39  throw cms::Exception("WrongConfig") << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
40  << "It has to be larger than 4 or can be set to -1 to take all jets.";
41 
42  produces<std::vector<std::vector<int>>>();
43  produces<int>("NumberOfConsideredJets");
44 }
45 
47  auto pOut = std::make_unique<std::vector<std::vector<int>>>();
48  auto pJetsConsidered = std::make_unique<int>(0);
49 
50  std::vector<int> match;
51  for (unsigned int i = 0; i < 4; ++i)
52  match.push_back(-1);
53 
54  // get jets
55  const auto& jets = evt.get(jetsToken_);
56 
57  // get leptons
58  const auto& leps = evt.get(lepsToken_);
59 
60  // skip events without lepton candidate or less than 4 jets or no MET
61  if (leps.empty() || jets.size() < 4) {
62  pOut->push_back(match);
63  evt.put(std::move(pOut));
64  *pJetsConsidered = jets.size();
65  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
66  return;
67  }
68 
69  unsigned maxNJets = maxNJets_;
70  if (maxNJets_ == -1 || (int)jets.size() < maxNJets_)
71  maxNJets = jets.size();
72  *pJetsConsidered = maxNJets;
73  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
74 
75  std::vector<bool> isBJet;
76  std::vector<bool> isLJet;
77  int cntBJets = 0;
78  if (useBTagging_) {
79  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
80  isBJet.push_back((jets[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_));
81  isLJet.push_back((jets[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_));
82  if (jets[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_)
83  cntBJets++;
84  }
85  }
86 
87  // -----------------------------------------------------
88  // associate those jets with maximum pt of the vectorial
89  // sum to the hadronic decay chain
90  // -----------------------------------------------------
91  double maxPt = -1.;
92  std::vector<int> maxPtIndices;
93  maxPtIndices.push_back(-1);
94  maxPtIndices.push_back(-1);
95  maxPtIndices.push_back(-1);
96  for (unsigned idx = 0; idx < maxNJets; ++idx) {
97  if (useBTagging_ && (!isLJet[idx] || (cntBJets <= 2 && isBJet[idx])))
98  continue;
99  for (unsigned jdx = (idx + 1); jdx < maxNJets; ++jdx) {
100  if (jdx == idx || (useBTagging_ && (!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) ||
101  (cntBJets == 3 && isBJet[idx] && isBJet[jdx]))))
102  continue;
103  for (unsigned kdx = 0; kdx < maxNJets; ++kdx) {
104  if (kdx == idx || kdx == jdx || (useBTagging_ && !isBJet[kdx]))
105  continue;
106  reco::Particle::LorentzVector sum = jets[idx].p4() + jets[jdx].p4() + jets[kdx].p4();
107  if (maxPt < 0. || maxPt < sum.pt()) {
108  maxPt = sum.pt();
109  maxPtIndices.clear();
110  maxPtIndices.push_back(idx);
111  maxPtIndices.push_back(jdx);
112  maxPtIndices.push_back(kdx);
113  }
114  }
115  }
116  }
117 
118  // -----------------------------------------------------
119  // associate those jets that get closest to the W mass
120  // with their invariant mass to the W boson
121  // -----------------------------------------------------
122  double wDist = -1.;
123  std::vector<int> closestToWMassIndices;
124  closestToWMassIndices.push_back(-1);
125  closestToWMassIndices.push_back(-1);
126  if (isValid(maxPtIndices[0], jets) && isValid(maxPtIndices[1], jets) && isValid(maxPtIndices[2], jets)) {
127  for (unsigned idx = 0; idx < maxPtIndices.size(); ++idx) {
128  for (unsigned jdx = 0; jdx < maxPtIndices.size(); ++jdx) {
129  if (jdx == idx || maxPtIndices[idx] > maxPtIndices[jdx] ||
130  (useBTagging_ &&
131  (!isLJet[maxPtIndices[idx]] || !isLJet[maxPtIndices[jdx]] ||
132  (cntBJets <= 2 && isBJet[maxPtIndices[idx]]) || (cntBJets <= 2 && isBJet[maxPtIndices[jdx]]) ||
133  (cntBJets == 3 && isBJet[maxPtIndices[idx]] && isBJet[maxPtIndices[jdx]]))))
134  continue;
135  reco::Particle::LorentzVector sum = jets[maxPtIndices[idx]].p4() + jets[maxPtIndices[jdx]].p4();
136  if (wDist < 0. || wDist > fabs(sum.mass() - wMass_)) {
137  wDist = fabs(sum.mass() - wMass_);
138  closestToWMassIndices.clear();
139  closestToWMassIndices.push_back(maxPtIndices[idx]);
140  closestToWMassIndices.push_back(maxPtIndices[jdx]);
141  }
142  }
143  }
144  }
145  int hadB = -1;
146  for (unsigned idx = 0; idx < maxPtIndices.size(); ++idx) {
147  // if this idx is not yet contained in the list of W mass candidates...
148  if (std::find(closestToWMassIndices.begin(), closestToWMassIndices.end(), maxPtIndices[idx]) ==
149  closestToWMassIndices.end()) {
150  hadB = maxPtIndices[idx];
151  break; // there should be no other cadidates!
152  }
153  }
154 
155  // -----------------------------------------------------
156  // associate the remaining jet with maximum pt of the
157  // vectorial sum with the leading lepton with the
158  // leptonic decay chain
159  // -----------------------------------------------------
160  maxPt = -1.;
161  int lepB = -1;
162  for (unsigned idx = 0; idx < maxNJets; ++idx) {
163  if (useBTagging_ && !isBJet[idx])
164  continue;
165  // make sure it's not used up already from the hadronic decay chain
166  if (std::find(maxPtIndices.begin(), maxPtIndices.end(), idx) == maxPtIndices.end()) {
167  reco::Particle::LorentzVector sum = jets[idx].p4() + leps[0].p4();
168  if (maxPt < 0. || maxPt < sum.pt()) {
169  maxPt = sum.pt();
170  lepB = idx;
171  }
172  }
173  }
174 
175  match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0];
176  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
179 
180  pOut->push_back(match);
181  evt.put(std::move(pOut));
182 }
183 
bool isValid(const int &idx, const std::vector< pat::Jet > &jets) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:347
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
Definition: HeavyIon.h:7
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
Definition: Jet.py:1
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
maxPt
Definition: PV_cfg.py:223
fixed size matrix
TtSemiLepJetCombMaxSumPtWMass(const edm::ParameterSet &)
HLT enums.
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &setup) const override
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
def move(src, dest)
Definition: eostools.py:511