CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepJetCombMaxSumPtWMass.cc
Go to the documentation of this file.
2 
4 
6 
8  jetsToken_ (consumes< std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets" ))),
9  lepsToken_ (consumes< edm::View<reco::RecoCandidate> >(cfg.getParameter<edm::InputTag>("leps" ))),
10  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
11  wMass_ (cfg.getParameter<double> ("wMass" )),
12  useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
13  bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )),
14  minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
15  maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets"))
16 {
17  if(maxNJets_<4 && maxNJets_!=-1)
18  throw cms::Exception("WrongConfig")
19  << "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 {
28 }
29 
30 void
32 {
33  std::auto_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
34  std::auto_ptr<int> pJetsConsidered(new 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.getByToken(jetsToken_, jets);
43 
44  // get leptons
46  evt.getByToken(lepsToken_, leps);
47 
48 
49  // skip events without lepton candidate or less than 4 jets or no MET
50  if(leps->empty() || jets->size() < 4){
51  pOut->push_back( match );
52  evt.put(pOut);
53  *pJetsConsidered = jets->size();
54  evt.put(pJetsConsidered, "NumberOfConsideredJets");
55  return;
56  }
57 
58  unsigned maxNJets = maxNJets_;
59  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
60  *pJetsConsidered = maxNJets;
61  evt.put(pJetsConsidered, "NumberOfConsideredJets");
62 
63  std::vector<bool> isBJet;
64  std::vector<bool> isLJet;
65  int cntBJets = 0;
66  if(useBTagging_) {
67  for(unsigned int idx=0; idx<maxNJets; ++idx) {
68  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
69  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
70  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
71  }
72  }
73 
74  // -----------------------------------------------------
75  // associate those jets with maximum pt of the vectorial
76  // sum to the hadronic decay chain
77  // -----------------------------------------------------
78  double maxPt=-1.;
79  std::vector<int> maxPtIndices;
80  maxPtIndices.push_back(-1);
81  maxPtIndices.push_back(-1);
82  maxPtIndices.push_back(-1);
83  for(unsigned idx=0; idx<maxNJets; ++idx){
84  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
85  for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){
86  if(jdx==idx || (useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx])))) continue;
87  for(unsigned kdx=0; kdx<maxNJets; ++kdx){
88  if(kdx==idx || kdx==jdx || (useBTagging_ && !isBJet[kdx])) continue;
90  (*jets)[idx].p4()+
91  (*jets)[jdx].p4()+
92  (*jets)[kdx].p4();
93  if( maxPt<0. || maxPt<sum.pt() ){
94  maxPt=sum.pt();
95  maxPtIndices.clear();
96  maxPtIndices.push_back(idx);
97  maxPtIndices.push_back(jdx);
98  maxPtIndices.push_back(kdx);
99  }
100  }
101  }
102  }
103 
104  // -----------------------------------------------------
105  // associate those jets that get closest to the W mass
106  // with their invariant mass to the W boson
107  // -----------------------------------------------------
108  double wDist =-1.;
109  std::vector<int> closestToWMassIndices;
110  closestToWMassIndices.push_back(-1);
111  closestToWMassIndices.push_back(-1);
112  if( isValid(maxPtIndices[0], jets) && isValid(maxPtIndices[1], jets) && isValid(maxPtIndices[2], jets)) {
113  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){
114  for(unsigned jdx=0; jdx<maxPtIndices.size(); ++jdx){
115  if( jdx==idx || maxPtIndices[idx]>maxPtIndices[jdx] || (useBTagging_ && (!isLJet[maxPtIndices[idx]] || !isLJet[maxPtIndices[jdx]] || (cntBJets<=2 && isBJet[maxPtIndices[idx]]) || (cntBJets<=2 && isBJet[maxPtIndices[jdx]]) || (cntBJets==3 && isBJet[maxPtIndices[idx]] && isBJet[maxPtIndices[jdx]])))) continue;
117  (*jets)[maxPtIndices[idx]].p4()+
118  (*jets)[maxPtIndices[jdx]].p4();
119  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
120  wDist=fabs(sum.mass()-wMass_);
121  closestToWMassIndices.clear();
122  closestToWMassIndices.push_back(maxPtIndices[idx]);
123  closestToWMassIndices.push_back(maxPtIndices[jdx]);
124  }
125  }
126  }
127  }
128  int hadB=-1;
129  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){
130  // if this idx is not yet contained in the list of W mass candidates...
131  if( std::find( closestToWMassIndices.begin(), closestToWMassIndices.end(), maxPtIndices[idx]) == closestToWMassIndices.end() ){
132  hadB = maxPtIndices[idx];
133  break; // there should be no other cadidates!
134  }
135  }
136 
137  // -----------------------------------------------------
138  // associate the remaining jet with maximum pt of the
139  // vectorial sum with the leading lepton with the
140  // leptonic decay chain
141  // -----------------------------------------------------
142  maxPt=-1.;
143  int lepB=-1;
144  for(unsigned idx=0; idx<maxNJets; ++idx){
145  if(useBTagging_ && !isBJet[idx]) continue;
146  // make sure it's not used up already from the hadronic decay chain
147  if( std::find(maxPtIndices.begin(), maxPtIndices.end(), idx) == maxPtIndices.end() ){
149  (*jets)[idx].p4()+(*leps)[ 0 ].p4();
150  if( maxPt<0. || maxPt<sum.pt() ){
151  maxPt=sum.pt();
152  lepB=idx;
153  }
154  }
155  }
156 
157  match[TtSemiLepEvtPartons::LightQ ] = closestToWMassIndices[0];
158  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
159  match[TtSemiLepEvtPartons::HadB ] = hadB;
160  match[TtSemiLepEvtPartons::LepB ] = lepB;
161 
162  pOut->push_back( match );
163  evt.put(pOut);
164 }
int i
Definition: DBlmapReader.cc:9
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
tuple cfg
Definition: looper.py:293
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
virtual void produce(edm::Event &evt, const edm::EventSetup &setup)
TtSemiLepJetCombMaxSumPtWMass(const edm::ParameterSet &)
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
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")