CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepJetCombWMassMaxSumPt.cc
Go to the documentation of this file.
2 
5 
7  jets_ (cfg.getParameter<edm::InputTag>("jets" )),
8  leps_ (cfg.getParameter<edm::InputTag>("leps" )),
9  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
10  wMass_ (cfg.getParameter<double> ("wMass" )),
11  useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
12  bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )),
13  minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
14  maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets"))
15 {
16  if(maxNJets_<4 && maxNJets_!=-1)
17  throw cms::Exception("WrongConfig")
18  << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
19  << "It has to be larger than 4 or can be set to -1 to take all jets.";
20 
21  produces<std::vector<std::vector<int> > >();
22  produces<int>("NumberOfConsideredJets");
23 }
24 
26 {
27 }
28 
29 void
31 {
32  std::auto_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
33  std::auto_ptr<int> pJetsConsidered(new int);
34 
35  std::vector<int> match;
36  for(unsigned int i = 0; i < 4; ++i)
37  match.push_back( -1 );
38 
39  // get jets
41  evt.getByLabel(jets_, jets);
42 
43  // get leptons
45  evt.getByLabel(leps_, leps);
46 
47 
48  // skip events without lepton candidate or less than 4 jets or no MET
49  if(leps->empty() || jets->size() < 4){
50  pOut->push_back( match );
51  evt.put(pOut);
52  *pJetsConsidered = jets->size();
53  evt.put(pJetsConsidered, "NumberOfConsideredJets");
54  return;
55  }
56 
57  unsigned maxNJets = maxNJets_;
58  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
59  *pJetsConsidered = maxNJets;
60  evt.put(pJetsConsidered, "NumberOfConsideredJets");
61 
62  std::vector<bool> isBJet;
63  std::vector<bool> isLJet;
64  int cntBJets = 0;
65  if(useBTagging_) {
66  for(unsigned int idx=0; idx<maxNJets; ++idx) {
67  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
68  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
69  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
70  }
71  }
72 
73  // -----------------------------------------------------
74  // associate those jets that get closest to the W mass
75  // with their invariant mass to the hadronic W boson
76  // -----------------------------------------------------
77  double wDist =-1.;
78  std::vector<int> closestToWMassIndices;
79  closestToWMassIndices.push_back(-1);
80  closestToWMassIndices.push_back(-1);
81  for(unsigned idx=0; idx<maxNJets; ++idx){
82  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
83  for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){
84  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
86  (*jets)[idx].p4()+
87  (*jets)[jdx].p4();
88  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
89  wDist=fabs(sum.mass()-wMass_);
90  closestToWMassIndices.clear();
91  closestToWMassIndices.push_back(idx);
92  closestToWMassIndices.push_back(jdx);
93  }
94  }
95  }
96 
97  // -----------------------------------------------------
98  // associate those jets with maximum pt of the vectorial
99  // sum to the hadronic decay chain
100  // -----------------------------------------------------
101  double maxPt=-1.;
102  int hadB=-1;
103  if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
104  for(unsigned idx=0; idx<maxNJets; ++idx){
105  if(useBTagging_ && !isBJet[idx]) continue;
106  // make sure it's not used up already from the hadronic W
107  if( (int)idx!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] ){
109  (*jets)[closestToWMassIndices[0]].p4()+
110  (*jets)[closestToWMassIndices[1]].p4()+
111  (*jets)[idx].p4();
112  if( maxPt<0. || maxPt<sum.pt() ){
113  maxPt=sum.pt();
114  hadB=idx;
115  }
116  }
117  }
118  }
119 
120  // -----------------------------------------------------
121  // associate the remaining jet with maximum pt of the
122  // vectorial sum with the leading lepton with the
123  // leptonic b quark
124  // -----------------------------------------------------
125  maxPt=-1.;
126  int lepB=-1;
127  for(unsigned idx=0; idx<maxNJets; ++idx){
128  if(useBTagging_ && !isBJet[idx]) continue;
129  // make sure it's not used up already from the hadronic decay chain
130  if( (int)idx!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] && (int)idx!=hadB) {
132  (*jets)[idx].p4()+(*leps)[ 0 ].p4();
133  if( maxPt<0. || maxPt<sum.pt() ){
134  maxPt=sum.pt();
135  lepB=idx;
136  }
137  }
138  }
139 
140  match[TtSemiLepEvtPartons::LightQ ] = closestToWMassIndices[0];
141  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
142  match[TtSemiLepEvtPartons::HadB ] = hadB;
143  match[TtSemiLepEvtPartons::LepB ] = lepB;
144 
145  pOut->push_back( match );
146  evt.put(pOut);
147 }
int i
Definition: DBlmapReader.cc:9
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
TtSemiLepJetCombWMassMaxSumPt(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
vector< PseudoJet > jets
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
virtual void produce(edm::Event &evt, const edm::EventSetup &setup)
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:25
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")