CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepJetCombGeom.cc
Go to the documentation of this file.
2 
3 #include "Math/VectorUtil.h"
4 
8 
10  jets_ (cfg.getParameter<edm::InputTag>("jets" )),
11  leps_ (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 }
26 
28 {
29 }
30 
31 void
33 {
34  std::auto_ptr< std::vector<std::vector<int> > >pOut (new std::vector<std::vector<int> >);
35 
36  std::vector<int> match;
37  for(unsigned int i = 0; i < 4; ++i)
38  match.push_back( -1 );
39 
40  // get jets
42  evt.getByLabel(jets_, jets);
43 
44  // get leptons
46  evt.getByLabel(leps_, leps);
47 
48  // skip events with without lepton candidate or less than 4 jets
49  if(leps->empty() || jets->size() < 4){
50  pOut->push_back( match );
51  evt.put(pOut);
52  return;
53  }
54 
55  unsigned maxNJets = maxNJets_;
56  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
57 
58  std::vector<bool> isBJet;
59  std::vector<bool> isLJet;
60  int cntBJets = 0;
61  if(useBTagging_) {
62  for(unsigned int idx=0; idx<maxNJets; ++idx) {
63  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
64  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
65  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
66  }
67  }
68 
69  // -----------------------------------------------------
70  // associate those two jets to the hadronic W boson that
71  // have the smallest distance to each other
72  // -----------------------------------------------------
73  double minDist=-1.;
74  int lightQ =-1;
75  int lightQBar=-1;
76  for(unsigned int idx=0; idx<maxNJets; ++idx){
77  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
78  for(unsigned int jdx=(idx+1); jdx<maxNJets; ++jdx){
79  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
80  double dist = distance((*jets)[idx].p4(), (*jets)[jdx].p4());
81  if( minDist<0. || dist<minDist ){
82  minDist=dist;
83  lightQ =idx;
84  lightQBar=jdx;
85  }
86  }
87  }
88 
90  if( isValid(lightQ, jets) && isValid(lightQBar, jets) )
91  wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4();
92 
93  // -----------------------------------------------------
94  // associate to the hadronic b quark the remaining jet
95  // that has the smallest distance to the hadronic W
96  // -----------------------------------------------------
97  minDist=-1.;
98  int hadB=-1;
99  if( isValid(lightQ, jets) && isValid(lightQBar, jets) ) {
100  for(unsigned int idx=0; idx<maxNJets; ++idx){
101  if(useBTagging_ && !isBJet[idx]) continue;
102  // make sure it's not used up already from the hadronic W
103  if( (int)idx!=lightQ && (int)idx!=lightQBar ) {
104  double dist = distance((*jets)[idx].p4(), wHad);
105  if( minDist<0. || dist<minDist ){
106  minDist=dist;
107  hadB=idx;
108  }
109  }
110  }
111  }
112 
113  // -----------------------------------------------------
114  // associate to the leptonic b quark the remaining jet
115  // that has the smallest distance to the leading lepton
116  // -----------------------------------------------------
117  minDist=-1.;
118  int lepB=-1;
119  for(unsigned int idx=0; idx<maxNJets; ++idx){
120  if(useBTagging_ && !isBJet[idx]) continue;
121  // make sure it's not used up already from the hadronic decay chain
122  if( (int)idx!=lightQ && (int)idx!=lightQBar && (int)idx!=hadB ){
123  double dist = distance((*jets)[idx].p4(), (*leps)[0].p4());
124  if( minDist<0. || dist<minDist ){
125  minDist=dist;
126  lepB=idx;
127  }
128  }
129  }
130 
131  match[TtSemiLepEvtPartons::LightQ ] = lightQ;
132  match[TtSemiLepEvtPartons::LightQBar] = lightQBar;
133  match[TtSemiLepEvtPartons::HadB ] = hadB;
134  match[TtSemiLepEvtPartons::LepB ] = lepB;
135 
136  pOut->push_back( match );
137  evt.put(pOut);
138 }
139 
140 double
142 {
143  // calculate the distance between two lorentz vectors
144  // using DeltaR or DeltaTheta
145  if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2);
146  return fabs(v1.theta() - v2.theta());
147 }
int i
Definition: DBlmapReader.cc:9
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
virtual void produce(edm::Event &evt, const edm::EventSetup &setup)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
TtSemiLepJetCombGeom(const edm::ParameterSet &)
double p4[4]
Definition: TauolaWrapper.h:92
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
double distance(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:26