CMS 3D CMS Logo

TtSemiLepJetCombGeom.cc
Go to the documentation of this file.
4 
9 
10 #include "Math/VectorUtil.h"
11 
13 public:
14  explicit TtSemiLepJetCombGeom(const edm::ParameterSet&);
15 
16 private:
17  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& setup) const override;
18 
19  bool isValid(const int idx, const std::vector<pat::Jet>& jets) const { return (0 <= idx && idx < (int)jets.size()); };
20  double distance(const math::XYZTLorentzVector&, const math::XYZTLorentzVector&) const;
21 
24  int maxNJets_;
25  bool useDeltaR_;
30 };
31 
33  : jetsToken_(consumes<std::vector<pat::Jet>>(cfg.getParameter<edm::InputTag>("jets"))),
34  lepsToken_(consumes<edm::View<reco::RecoCandidate>>(cfg.getParameter<edm::InputTag>("leps"))),
35  maxNJets_(cfg.getParameter<int>("maxNJets")),
36  useDeltaR_(cfg.getParameter<bool>("useDeltaR")),
37  useBTagging_(cfg.getParameter<bool>("useBTagging")),
38  bTagAlgorithm_(cfg.getParameter<std::string>("bTagAlgorithm")),
39  minBDiscBJets_(cfg.getParameter<double>("minBDiscBJets")),
40  maxBDiscLightJets_(cfg.getParameter<double>("maxBDiscLightJets")) {
41  if (maxNJets_ < 4 && maxNJets_ != -1)
42  throw cms::Exception("WrongConfig") << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
43  << "It has to be larger than 4 or can be set to -1 to take all jets.";
44 
45  produces<std::vector<std::vector<int>>>();
46  produces<int>("NumberOfConsideredJets");
47 }
48 
50  auto pOut = std::make_unique<std::vector<std::vector<int>>>();
51  auto pJetsConsidered = std::make_unique<int>(0);
52 
53  std::vector<int> match;
54  for (unsigned int i = 0; i < 4; ++i)
55  match.push_back(-1);
56 
57  // get jets
58  const auto& jets = evt.get(jetsToken_);
59 
60  // get leptons
61  const auto& leps = evt.get(lepsToken_);
62 
63  // skip events without lepton candidate or less than 4 jets
64  if (leps.empty() || jets.size() < 4) {
65  pOut->push_back(match);
66  evt.put(std::move(pOut));
67  *pJetsConsidered = jets.size();
68  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
69  return;
70  }
71 
72  unsigned maxNJets = maxNJets_;
73  if (maxNJets_ == -1 || (int)jets.size() < maxNJets_)
74  maxNJets = jets.size();
75  *pJetsConsidered = maxNJets;
76  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
77 
78  std::vector<bool> isBJet;
79  std::vector<bool> isLJet;
80  int cntBJets = 0;
81  if (useBTagging_) {
82  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
83  isBJet.push_back((jets[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_));
84  isLJet.push_back((jets[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_));
85  if (jets[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_)
86  cntBJets++;
87  }
88  }
89 
90  // -----------------------------------------------------
91  // associate those two jets to the hadronic W boson that
92  // have the smallest distance to each other
93  // -----------------------------------------------------
94  double minDist = -1.;
95  int lightQ = -1;
96  int lightQBar = -1;
97  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
98  if (useBTagging_ && (!isLJet[idx] || (cntBJets <= 2 && isBJet[idx])))
99  continue;
100  for (unsigned int jdx = (idx + 1); jdx < maxNJets; ++jdx) {
101  if (useBTagging_ &&
102  (!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) || (cntBJets == 3 && isBJet[idx] && isBJet[jdx])))
103  continue;
104  double dist = distance(jets[idx].p4(), jets[jdx].p4());
105  if (minDist < 0. || dist < minDist) {
106  minDist = dist;
107  lightQ = idx;
108  lightQBar = jdx;
109  }
110  }
111  }
112 
114  if (isValid(lightQ, jets) && isValid(lightQBar, jets))
115  wHad = jets[lightQ].p4() + jets[lightQBar].p4();
116 
117  // -----------------------------------------------------
118  // associate to the hadronic b quark the remaining jet
119  // that has the smallest distance to the hadronic W
120  // -----------------------------------------------------
121  minDist = -1.;
122  int hadB = -1;
123  if (isValid(lightQ, jets) && isValid(lightQBar, jets)) {
124  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
125  if (useBTagging_ && !isBJet[idx])
126  continue;
127  // make sure it's not used up already from the hadronic W
128  if ((int)idx != lightQ && (int)idx != lightQBar) {
129  double dist = distance(jets[idx].p4(), wHad);
130  if (minDist < 0. || dist < minDist) {
131  minDist = dist;
132  hadB = idx;
133  }
134  }
135  }
136  }
137 
138  // -----------------------------------------------------
139  // associate to the leptonic b quark the remaining jet
140  // that has the smallest distance to the leading lepton
141  // -----------------------------------------------------
142  minDist = -1.;
143  int lepB = -1;
144  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
145  if (useBTagging_ && !isBJet[idx])
146  continue;
147  // make sure it's not used up already from the hadronic decay chain
148  if ((int)idx != lightQ && (int)idx != lightQBar && (int)idx != hadB) {
149  double dist = distance(jets[idx].p4(), leps[0].p4());
150  if (minDist < 0. || dist < minDist) {
151  minDist = dist;
152  lepB = idx;
153  }
154  }
155  }
156 
161 
162  pOut->push_back(match);
163  evt.put(std::move(pOut));
164 }
165 
167  // calculate the distance between two lorentz vectors
168  // using DeltaR or DeltaTheta
169  if (useDeltaR_)
170  return ROOT::Math::VectorUtil::DeltaR(v1, v2);
171  return fabs(v1.theta() - v2.theta());
172 }
173 
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:346
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double distance(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &) const
Definition: HeavyIon.h:7
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool isValid(const int idx, const std::vector< pat::Jet > &jets) const
Definition: Jet.py:1
TtSemiLepJetCombGeom(const edm::ParameterSet &)
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &setup) const override
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
fixed size matrix
HLT enums.
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