CMS 3D CMS Logo

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 {
19  if(maxNJets_<4 && maxNJets_!=-1)
20  throw cms::Exception("WrongConfig")
21  << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
22  << "It has to be larger than 4 or can be set to -1 to take all jets.";
23 
24  produces<std::vector<std::vector<int> > >();
25  produces<int>("NumberOfConsideredJets");
26 }
27 
29 {
30 }
31 
32 void
34 {
35  std::unique_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
36  std::unique_ptr<int> pJetsConsidered(new int);
37 
38  std::vector<int> match;
39  for(unsigned int i = 0; i < 4; ++i)
40  match.push_back( -1 );
41 
42  // get jets
44  evt.getByToken(jetsToken_, jets);
45 
46  // get leptons
48  evt.getByToken(lepsToken_, leps);
49 
50  // skip events without lepton candidate or less than 4 jets
51  if(leps->empty() || jets->size() < 4){
52  pOut->push_back( match );
53  evt.put(std::move(pOut));
54  *pJetsConsidered = jets->size();
55  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
56  return;
57  }
58 
59  unsigned maxNJets = maxNJets_;
60  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
61  *pJetsConsidered = maxNJets;
62  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
63 
64  std::vector<bool> isBJet;
65  std::vector<bool> isLJet;
66  int cntBJets = 0;
67  if(useBTagging_) {
68  for(unsigned int idx=0; idx<maxNJets; ++idx) {
69  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
70  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
71  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
72  }
73  }
74 
75  // -----------------------------------------------------
76  // associate those two jets to the hadronic W boson that
77  // have the smallest distance to each other
78  // -----------------------------------------------------
79  double minDist=-1.;
80  int lightQ =-1;
81  int lightQBar=-1;
82  for(unsigned int idx=0; idx<maxNJets; ++idx){
83  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
84  for(unsigned int jdx=(idx+1); jdx<maxNJets; ++jdx){
85  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
86  double dist = distance((*jets)[idx].p4(), (*jets)[jdx].p4());
87  if( minDist<0. || dist<minDist ){
88  minDist=dist;
89  lightQ =idx;
90  lightQBar=jdx;
91  }
92  }
93  }
94 
96  if( isValid(lightQ, jets) && isValid(lightQBar, jets) )
97  wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4();
98 
99  // -----------------------------------------------------
100  // associate to the hadronic b quark the remaining jet
101  // that has the smallest distance to the hadronic W
102  // -----------------------------------------------------
103  minDist=-1.;
104  int hadB=-1;
105  if( isValid(lightQ, jets) && isValid(lightQBar, jets) ) {
106  for(unsigned int idx=0; idx<maxNJets; ++idx){
107  if(useBTagging_ && !isBJet[idx]) 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]) continue;
127  // make sure it's not used up already from the hadronic decay chain
128  if( (int)idx!=lightQ && (int)idx!=lightQBar && (int)idx!=hadB ){
129  double dist = distance((*jets)[idx].p4(), (*leps)[0].p4());
130  if( minDist<0. || dist<minDist ){
131  minDist=dist;
132  lepB=idx;
133  }
134  }
135  }
136 
137  match[TtSemiLepEvtPartons::LightQ ] = lightQ;
138  match[TtSemiLepEvtPartons::LightQBar] = lightQBar;
139  match[TtSemiLepEvtPartons::HadB ] = hadB;
140  match[TtSemiLepEvtPartons::LepB ] = lepB;
141 
142  pOut->push_back( match );
143  evt.put(std::move(pOut));
144 }
145 
146 double
148 {
149  // calculate the distance between two lorentz vectors
150  // using DeltaR or DeltaTheta
151  if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2);
152  return fabs(v1.theta() - v2.theta());
153 }
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: HeavyIon.h:7
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Definition: Jet.py:1
TtSemiLepJetCombGeom(const edm::ParameterSet &)
double p4[4]
Definition: TauolaWrapper.h:92
vector< PseudoJet > jets
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
maxNJets
maximum number of jets taken into account per event for each hypothesis (this parameter is used in th...
double distance(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
fixed size matrix
HLT enums.
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
def move(src, dest)
Definition: eostools.py:510