CMS 3D CMS Logo

TtSemiLepJetCombWMassMaxSumPt.cc
Go to the documentation of this file.
2 
4 
6  jetsToken_ (consumes< std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets" ))),
7  lepsToken_ (consumes< edm::View<reco::RecoCandidate> >(cfg.getParameter<edm::InputTag>("leps" ))),
8  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
9  wMass_ (cfg.getParameter<double> ("wMass" )),
10  useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
11  bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )),
12  minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
13  maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets"))
14 {
15  if(maxNJets_<4 && maxNJets_!=-1)
16  throw cms::Exception("WrongConfig")
17  << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
18  << "It has to be larger than 4 or can be set to -1 to take all jets.";
19 
20  produces<std::vector<std::vector<int> > >();
21  produces<int>("NumberOfConsideredJets");
22 }
23 
25 {
26 }
27 
28 void
30 {
31  std::unique_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
32  std::unique_ptr<int> pJetsConsidered(new int);
33 
34  std::vector<int> match;
35  for(unsigned int i = 0; i < 4; ++i)
36  match.push_back( -1 );
37 
38  // get jets
40  evt.getByToken(jetsToken_, jets);
41 
42  // get leptons
44  evt.getByToken(lepsToken_, leps);
45 
46 
47  // skip events without lepton candidate or less than 4 jets or no MET
48  if(leps->empty() || jets->size() < 4){
49  pOut->push_back( match );
50  evt.put(std::move(pOut));
51  *pJetsConsidered = jets->size();
52  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
53  return;
54  }
55 
56  unsigned maxNJets = maxNJets_;
57  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
58  *pJetsConsidered = maxNJets;
59  evt.put(std::move(pJetsConsidered), "NumberOfConsideredJets");
60 
61  std::vector<bool> isBJet;
62  std::vector<bool> isLJet;
63  int cntBJets = 0;
64  if(useBTagging_) {
65  for(unsigned int idx=0; idx<maxNJets; ++idx) {
66  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
67  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
68  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
69  }
70  }
71 
72  // -----------------------------------------------------
73  // associate those jets that get closest to the W mass
74  // with their invariant mass to the hadronic W boson
75  // -----------------------------------------------------
76  double wDist =-1.;
77  std::vector<int> closestToWMassIndices;
78  closestToWMassIndices.push_back(-1);
79  closestToWMassIndices.push_back(-1);
80  for(unsigned idx=0; idx<maxNJets; ++idx){
81  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
82  for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){
83  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
85  (*jets)[idx].p4()+
86  (*jets)[jdx].p4();
87  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
88  wDist=fabs(sum.mass()-wMass_);
89  closestToWMassIndices.clear();
90  closestToWMassIndices.push_back(idx);
91  closestToWMassIndices.push_back(jdx);
92  }
93  }
94  }
95 
96  // -----------------------------------------------------
97  // associate those jets with maximum pt of the vectorial
98  // sum to the hadronic decay chain
99  // -----------------------------------------------------
100  double maxPt=-1.;
101  int hadB=-1;
102  if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
103  for(unsigned idx=0; idx<maxNJets; ++idx){
104  if(useBTagging_ && !isBJet[idx]) continue;
105  // make sure it's not used up already from the hadronic W
106  if( (int)idx!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] ){
108  (*jets)[closestToWMassIndices[0]].p4()+
109  (*jets)[closestToWMassIndices[1]].p4()+
110  (*jets)[idx].p4();
111  if( maxPt<0. || maxPt<sum.pt() ){
112  maxPt=sum.pt();
113  hadB=idx;
114  }
115  }
116  }
117  }
118 
119  // -----------------------------------------------------
120  // associate the remaining jet with maximum pt of the
121  // vectorial sum with the leading lepton with the
122  // leptonic b quark
123  // -----------------------------------------------------
124  maxPt=-1.;
125  int lepB=-1;
126  for(unsigned idx=0; idx<maxNJets; ++idx){
127  if(useBTagging_ && !isBJet[idx]) continue;
128  // make sure it's not used up already from the hadronic decay chain
129  if( (int)idx!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] && (int)idx!=hadB) {
131  (*jets)[idx].p4()+(*leps)[ 0 ].p4();
132  if( maxPt<0. || maxPt<sum.pt() ){
133  maxPt=sum.pt();
134  lepB=idx;
135  }
136  }
137  }
138 
139  match[TtSemiLepEvtPartons::LightQ ] = closestToWMassIndices[0];
140  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
141  match[TtSemiLepEvtPartons::HadB ] = hadB;
142  match[TtSemiLepEvtPartons::LepB ] = lepB;
143 
144  pOut->push_back( match );
145  evt.put(std::move(pOut));
146 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
TtSemiLepJetCombWMassMaxSumPt(const edm::ParameterSet &)
void produce(edm::Event &evt, const edm::EventSetup &setup) override
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: HeavyIon.h:7
Definition: Jet.py:1
vector< PseudoJet > jets
maxNJets
maximum number of jets taken into account per event for each hypothesis (this parameter is used in th...
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:510