CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepJetCombWMassDeltaTopMass.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  metsToken_ (consumes< std::vector<pat::MET> >(cfg.getParameter<edm::InputTag>("mets" ))),
11  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
12  wMass_ (cfg.getParameter<double> ("wMass" )),
13  useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
14  bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )),
15  minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
16  maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets")),
17  neutrinoSolutionType_(cfg.getParameter<int> ("neutrinoSolutionType"))
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::auto_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
36  std::auto_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  // get MET
52  evt.getByToken(metsToken_, mets);
53 
54  // skip events without lepton candidate or less than 4 jets or no MET
55  if(leps->empty() || jets->size() < 4 || mets->empty()){
56  pOut->push_back( match );
57  evt.put(pOut);
58  *pJetsConsidered = jets->size();
59  evt.put(pJetsConsidered, "NumberOfConsideredJets");
60  return;
61  }
62 
63  unsigned maxNJets = maxNJets_;
64  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
65  *pJetsConsidered = maxNJets;
66  evt.put(pJetsConsidered, "NumberOfConsideredJets");
67 
68  std::vector<bool> isBJet;
69  std::vector<bool> isLJet;
70  int cntBJets = 0;
71  if(useBTagging_) {
72  for(unsigned int idx=0; idx<maxNJets; ++idx) {
73  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
74  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
75  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
76  }
77  }
78 
79  // -----------------------------------------------------
80  // associate those jets that get closest to the W mass
81  // with their invariant mass to the hadronic W boson
82  // -----------------------------------------------------
83  double wDist =-1.;
84  std::vector<int> closestToWMassIndices;
85  closestToWMassIndices.push_back(-1);
86  closestToWMassIndices.push_back(-1);
87  for(unsigned idx=0; idx<maxNJets; ++idx){
88  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
89  for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){
90  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
92  (*jets)[idx].p4()+
93  (*jets)[jdx].p4();
94  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
95  wDist=fabs(sum.mass()-wMass_);
96  closestToWMassIndices.clear();
97  closestToWMassIndices.push_back(idx);
98  closestToWMassIndices.push_back(jdx);
99  }
100  }
101  }
102 
103  // -----------------------------------------------------
104  // build a leptonic W boson from the lepton and the MET
105  // -----------------------------------------------------
106  double neutrino_pz = 0.;
107  if(neutrinoSolutionType_!=-1) {
108  MEzCalculator mez;
109  mez.SetMET( *mets->begin() );
110  if( dynamic_cast<const reco::Muon*>(&(leps->front())) )
111  mez.SetLepton( (*leps)[0], true );
112  else if( dynamic_cast<const reco::GsfElectron*>(&(leps->front())) )
113  mez.SetLepton( (*leps)[0], false );
114  else
115  throw cms::Exception("UnimplementedFeature") << "Type of lepton given together with MET for solving neutrino kinematics is neither muon nor electron.\n";
116  neutrino_pz = mez.Calculate(neutrinoSolutionType_);
117  }
118  const math::XYZTLorentzVector neutrino( mets->begin()->px(), mets->begin()->py(), neutrino_pz, sqrt(mets->begin()->px()*mets->begin()->px()
119  + mets->begin()->py()*mets->begin()->py()
120  + neutrino_pz*neutrino_pz) );
121  const reco::Particle::LorentzVector lepW = neutrino + leps->front().p4();
122 
123  // -----------------------------------------------------
124  // associate those jets to the hadronic and the leptonic
125  // b quark that minimize the difference between
126  // hadronic and leptonic top mass
127  // -----------------------------------------------------
128  double deltaTop=-1.;
129  int hadB=-1;
130  int lepB=-1;
131  if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
132  const reco::Particle::LorentzVector hadW =
133  (*jets)[closestToWMassIndices[0]].p4()+
134  (*jets)[closestToWMassIndices[1]].p4();
135  // find hadronic b candidate
136  for(unsigned idx=0; idx<maxNJets; ++idx){
137  if(useBTagging_ && !isBJet[idx]) continue;
138  // make sure it's not used up already from the hadronic W
139  if( (int)idx!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] ){
140  reco::Particle::LorentzVector hadTop = hadW + (*jets)[idx].p4();
141  // find leptonic b candidate
142  for(unsigned jdx=0; jdx<maxNJets; ++jdx){
143  if(useBTagging_ && !isBJet[jdx]) continue;
144  // make sure it's not used up already from the hadronic branch
145  if( (int)jdx!=closestToWMassIndices[0] && (int)jdx!=closestToWMassIndices[1] && jdx!=idx ){
146  reco::Particle::LorentzVector lepTop = lepW + (*jets)[jdx].p4();
147  if( deltaTop<0. || deltaTop>fabs(hadTop.mass()-lepTop.mass()) ){
148  deltaTop=fabs(hadTop.mass()-lepTop.mass());
149  hadB=idx;
150  lepB=jdx;
151  }
152  }
153  }
154  }
155  }
156  }
157 
158  match[TtSemiLepEvtPartons::LightQ ] = closestToWMassIndices[0];
159  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
160  match[TtSemiLepEvtPartons::HadB ] = hadB;
161  match[TtSemiLepEvtPartons::LepB ] = lepB;
162 
163  pOut->push_back( match );
164  evt.put(pOut);
165 }
void SetLepton(const pat::Particle &lepton, bool isMuon=true)
Set lepton.
Definition: MEzCalculator.h:32
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:259
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
virtual void produce(edm::Event &evt, const edm::EventSetup &setup)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
void SetMET(const pat::MET &MET)
Set MET.
Definition: MEzCalculator.h:26
double Calculate(int type=1)
member functions
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
TtSemiLepJetCombWMassDeltaTopMass(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
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="")