CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtSemiLepJetCombGeom.cc
Go to the documentation of this file.
2 
3 #include "Math/VectorUtil.h"
4 
8 
10  : jetsToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
11  lepsToken_(consumes<edm::View<reco::RecoCandidate> >(cfg.getParameter<edm::InputTag>("leps"))),
12  maxNJets_(cfg.getParameter<int>("maxNJets")),
13  useDeltaR_(cfg.getParameter<bool>("useDeltaR")),
14  useBTagging_(cfg.getParameter<bool>("useBTagging")),
15  bTagAlgorithm_(cfg.getParameter<std::string>("bTagAlgorithm")),
16  minBDiscBJets_(cfg.getParameter<double>("minBDiscBJets")),
17  maxBDiscLightJets_(cfg.getParameter<double>("maxBDiscLightJets")) {
18  if (maxNJets_ < 4 && maxNJets_ != -1)
19  throw cms::Exception("WrongConfig") << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
20  << "It has to be larger than 4 or can be set to -1 to take all jets.";
21 
22  produces<std::vector<std::vector<int> > >();
23  produces<int>("NumberOfConsideredJets");
24 }
25 
27 
29  std::unique_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
30  std::unique_ptr<int> pJetsConsidered(new int);
31 
32  std::vector<int> match;
33  for (unsigned int i = 0; i < 4; ++i)
34  match.push_back(-1);
35 
36  // get jets
38  evt.getByToken(jetsToken_, jets);
39 
40  // get leptons
42  evt.getByToken(lepsToken_, leps);
43 
44  // skip events without lepton candidate or less than 4 jets
45  if (leps->empty() || jets->size() < 4) {
46  pOut->push_back(match);
47  evt.put(std::move(pOut));
48  *pJetsConsidered = jets->size();
49  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
50  return;
51  }
52 
53  unsigned maxNJets = maxNJets_;
54  if (maxNJets_ == -1 || (int)jets->size() < maxNJets_)
55  maxNJets = jets->size();
56  *pJetsConsidered = maxNJets;
57  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
58 
59  std::vector<bool> isBJet;
60  std::vector<bool> isLJet;
61  int cntBJets = 0;
62  if (useBTagging_) {
63  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
64  isBJet.push_back(((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_));
65  isLJet.push_back(((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_));
66  if ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_)
67  cntBJets++;
68  }
69  }
70 
71  // -----------------------------------------------------
72  // associate those two jets to the hadronic W boson that
73  // have the smallest distance to each other
74  // -----------------------------------------------------
75  double minDist = -1.;
76  int lightQ = -1;
77  int lightQBar = -1;
78  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
79  if (useBTagging_ && (!isLJet[idx] || (cntBJets <= 2 && isBJet[idx])))
80  continue;
81  for (unsigned int jdx = (idx + 1); jdx < maxNJets; ++jdx) {
82  if (useBTagging_ &&
83  (!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) || (cntBJets == 3 && isBJet[idx] && isBJet[jdx])))
84  continue;
85  double dist = distance((*jets)[idx].p4(), (*jets)[jdx].p4());
86  if (minDist < 0. || dist < minDist) {
87  minDist = dist;
88  lightQ = idx;
89  lightQBar = jdx;
90  }
91  }
92  }
93 
95  if (isValid(lightQ, jets) && isValid(lightQBar, jets))
96  wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4();
97 
98  // -----------------------------------------------------
99  // associate to the hadronic b quark the remaining jet
100  // that has the smallest distance to the hadronic W
101  // -----------------------------------------------------
102  minDist = -1.;
103  int hadB = -1;
104  if (isValid(lightQ, jets) && isValid(lightQBar, jets)) {
105  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
106  if (useBTagging_ && !isBJet[idx])
107  continue;
108  // make sure it's not used up already from the hadronic W
109  if ((int)idx != lightQ && (int)idx != lightQBar) {
110  double dist = distance((*jets)[idx].p4(), wHad);
111  if (minDist < 0. || dist < minDist) {
112  minDist = dist;
113  hadB = idx;
114  }
115  }
116  }
117  }
118 
119  // -----------------------------------------------------
120  // associate to the leptonic b quark the remaining jet
121  // that has the smallest distance to the leading lepton
122  // -----------------------------------------------------
123  minDist = -1.;
124  int lepB = -1;
125  for (unsigned int idx = 0; idx < maxNJets; ++idx) {
126  if (useBTagging_ && !isBJet[idx])
127  continue;
128  // make sure it's not used up already from the hadronic decay chain
129  if ((int)idx != lightQ && (int)idx != lightQBar && (int)idx != hadB) {
130  double dist = distance((*jets)[idx].p4(), (*leps)[0].p4());
131  if (minDist < 0. || dist < minDist) {
132  minDist = dist;
133  lepB = idx;
134  }
135  }
136  }
137 
138  match[TtSemiLepEvtPartons::LightQ] = lightQ;
139  match[TtSemiLepEvtPartons::LightQBar] = lightQBar;
140  match[TtSemiLepEvtPartons::HadB] = hadB;
141  match[TtSemiLepEvtPartons::LepB] = lepB;
142 
143  pOut->push_back(match);
144  evt.put(std::move(pOut));
145 }
146 
148  // calculate the distance between two lorentz vectors
149  // using DeltaR or DeltaTheta
150  if (useDeltaR_)
151  return ROOT::Math::VectorUtil::DeltaR(v1, v2);
152  return fabs(v1.theta() - v2.theta());
153 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple cfg
Definition: looper.py:296
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
TtSemiLepJetCombGeom(const edm::ParameterSet &)
vector< PseudoJet > jets
def move
Definition: eostools.py:511
double distance(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
constexpr char Jet[]
Definition: modules.cc:9
void produce(edm::Event &evt, const edm::EventSetup &setup) 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