CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepHypWMassMaxSumPt.cc
Go to the documentation of this file.
3 
5  TtSemiLepHypothesis( cfg ),
6  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
7  wMass_ (cfg.getParameter<double> ("wMass" )),
8  useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
9  bTagAlgorithm_ (cfg.getParameter<std::string>("bTagAlgorithm" )),
10  minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
11  maxBDiscLightJets_ (cfg.getParameter<double> ("maxBDiscLightJets" )),
12  neutrinoSolutionType_(cfg.getParameter<int> ("neutrinoSolutionType"))
13 {
14  if(maxNJets_<4 && maxNJets_!=-1)
15  throw cms::Exception("WrongConfig")
16  << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
17  << "It has to be larger than 4 or can be set to -1 to take all jets.";
18 }
19 
21 
22 void
25  const edm::Handle<std::vector<pat::MET> >& mets,
26  const edm::Handle<std::vector<pat::Jet> >& jets,
27  std::vector<int>& match, const unsigned int iComb)
28 {
29  if(leps->empty() || mets->empty() || jets->size()<4){
30  // create empty hypothesis
31  return;
32  }
33 
34  int maxNJets = maxNJets_;
35  if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
36 
37  std::vector<bool> isBJet;
38  std::vector<bool> isLJet;
39  int cntBJets = 0;
40  if(useBTagging_) {
41  for(int idx=0; idx<maxNJets; ++idx) {
42  isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
43  isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
44  if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
45  }
46  }
47 
48  match.clear();
49  for(unsigned int i=0; i<5; ++i)
50  match.push_back(-1);
51 
52  // -----------------------------------------------------
53  // associate those jets that get closest to the W mass
54  // with their invariant mass to the hadronic W boson
55  // -----------------------------------------------------
56  double wDist =-1.;
57  std::vector<int> closestToWMassIndices;
58  closestToWMassIndices.push_back(-1);
59  closestToWMassIndices.push_back(-1);
60  for(int idx=0; idx<maxNJets; ++idx){
61  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
62  for(int jdx=(idx+1); jdx<maxNJets; ++jdx){
63  if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
65  (*jets)[idx].p4()+
66  (*jets)[jdx].p4();
67  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
68  wDist=fabs(sum.mass()-wMass_);
69  closestToWMassIndices.clear();
70  closestToWMassIndices.push_back(idx);
71  closestToWMassIndices.push_back(jdx);
72  }
73  }
74  }
75 
76  // -----------------------------------------------------
77  // associate those jets with maximum pt of the vectorial
78  // sum to the hadronic decay chain
79  // -----------------------------------------------------
80  double maxPt=-1.;
81  int hadB=-1;
82  if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
83  for(int idx=0; idx<maxNJets; ++idx){
84  if(useBTagging_ && !isBJet[idx]) continue;
85  // make sure it's not used up already from the hadronic W
86  if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] ){
88  (*jets)[closestToWMassIndices[0]].p4()+
89  (*jets)[closestToWMassIndices[1]].p4()+
90  (*jets)[idx].p4();
91  if( maxPt<0. || maxPt<sum.pt() ){
92  maxPt=sum.pt();
93  hadB=idx;
94  }
95  }
96  }
97  }
98 
99  // -----------------------------------------------------
100  // associate the remaining jet with maximum pt of the
101  // vectorial sum with the leading lepton with the
102  // leptonic b quark
103  // -----------------------------------------------------
104  maxPt=-1.;
105  int lepB=-1;
106  for(int idx=0; idx<maxNJets; ++idx){
107  if(useBTagging_ && !isBJet[idx]) continue;
108  // make sure it's not used up already from the hadronic decay chain
109  if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] && idx!=hadB) {
111  (*jets)[idx].p4()+(*leps)[ 0 ].p4();
112  if( maxPt<0. || maxPt<sum.pt() ){
113  maxPt=sum.pt();
114  lepB=idx;
115  }
116  }
117  }
118 
119  // -----------------------------------------------------
120  // add jets
121  // -----------------------------------------------------
122  if( isValid(closestToWMassIndices[0], jets) ){
123  setCandidate(jets, closestToWMassIndices[0], lightQ_, jetCorrectionLevel("wQuarkMix"));
124  match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0];
125  }
126 
127  if( isValid(closestToWMassIndices[1], jets) ){
128  setCandidate(jets, closestToWMassIndices[1], lightQBar_, jetCorrectionLevel("wQuarkMix"));
129  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
130  }
131 
132  if( isValid(hadB, jets) ){
133  setCandidate(jets, hadB, hadronicB_, jetCorrectionLevel("bQuark"));
134  match[TtSemiLepEvtPartons::HadB] = hadB;
135  }
136 
137  if( isValid(lepB, jets) ){
138  setCandidate(jets, lepB, leptonicB_, jetCorrectionLevel("bQuark"));
139  match[TtSemiLepEvtPartons::LepB] = lepB;
140  }
141 
142  // -----------------------------------------------------
143  // add lepton
144  // -----------------------------------------------------
145  setCandidate(leps, 0, lepton_);
146  match[TtSemiLepEvtPartons::Lepton] = 0;
147 
148  // -----------------------------------------------------
149  // add neutrino
150  // -----------------------------------------------------
151  if(neutrinoSolutionType_ == -1)
152  setCandidate(mets, 0, neutrino_);
153  else
154  setNeutrino(mets, leps, 0, neutrinoSolutionType_);
155 }
int i
Definition: DBlmapReader.cc:9
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
reco::ShallowClonePtrCandidate * lepton_
reco::ShallowClonePtrCandidate * lightQBar_
reco::ShallowClonePtrCandidate * neutrino_
reco::ShallowClonePtrCandidate * hadronicB_
TtSemiLepHypWMassMaxSumPt(const edm::ParameterSet &)
void setNeutrino(const edm::Handle< std::vector< pat::MET > > &met, const edm::Handle< edm::View< reco::RecoCandidate > > &leps, const int &idx, const int &type)
set neutrino, using mW = 80.4 to calculate the neutrino pz
virtual void buildHypo(edm::Event &, const edm::Handle< edm::View< reco::RecoCandidate > > &, const edm::Handle< std::vector< pat::MET > > &, const edm::Handle< std::vector< pat::Jet > > &, std::vector< int > &, const unsigned int iComb)
build event hypothesis from the reco objects of a semi-leptonic event
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
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
reco::ShallowClonePtrCandidate * lightQ_
reco::ShallowClonePtrCandidate * leptonicB_
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...