CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepHypMaxSumPtWMass.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(int i=0; i<5; ++i)
50  match.push_back(-1);
51 
52  // -----------------------------------------------------
53  // associate those jets with maximum pt of the vectorial
54  // sum to the hadronic decay chain
55  // -----------------------------------------------------
56  double maxPt=-1.;
57  std::vector<int> maxPtIndices;
58  maxPtIndices.push_back(-1);
59  maxPtIndices.push_back(-1);
60  maxPtIndices.push_back(-1);
61  for(int idx=0; idx<maxNJets; ++idx){
62  if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
63  for(int jdx=(idx+1); jdx<maxNJets; ++jdx){
64  if(jdx==idx || (useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx])))) continue;
65  for(int kdx=0; kdx<maxNJets; ++kdx){
66  if(kdx==idx || kdx==jdx || (useBTagging_ && !isBJet[kdx])) continue;
68  (*jets)[idx].p4()+
69  (*jets)[jdx].p4()+
70  (*jets)[kdx].p4();
71  if( maxPt<0. || maxPt<sum.pt() ){
72  maxPt=sum.pt();
73  maxPtIndices.clear();
74  maxPtIndices.push_back(idx);
75  maxPtIndices.push_back(jdx);
76  maxPtIndices.push_back(kdx);
77  }
78  }
79  }
80  }
81 
82  // -----------------------------------------------------
83  // associate those jets that get closest to the W mass
84  // with their invariant mass to the W boson
85  // -----------------------------------------------------
86  double wDist =-1.;
87  std::vector<int> closestToWMassIndices;
88  closestToWMassIndices.push_back(-1);
89  closestToWMassIndices.push_back(-1);
90  if( isValid(maxPtIndices[0], jets) && isValid(maxPtIndices[1], jets) && isValid(maxPtIndices[2], jets)) {
91  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){
92  for(unsigned jdx=0; jdx<maxPtIndices.size(); ++jdx){
93  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;
95  (*jets)[maxPtIndices[idx]].p4()+
96  (*jets)[maxPtIndices[jdx]].p4();
97  if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
98  wDist=fabs(sum.mass()-wMass_);
99  closestToWMassIndices.clear();
100  closestToWMassIndices.push_back(maxPtIndices[idx]);
101  closestToWMassIndices.push_back(maxPtIndices[jdx]);
102  }
103  }
104  }
105  }
106 
107  // -----------------------------------------------------
108  // associate the remaining jet with maximum pt of the
109  // vectorial sum with the leading lepton with the
110  // leptonic decay chain
111  // -----------------------------------------------------
112  maxPt=-1.;
113  int lepB=-1;
114  for(int idx=0; idx<maxNJets; ++idx){
115  if(useBTagging_ && !isBJet[idx]) continue;
116  // make sure it's not used up already from the hadronic decay chain
117  if( std::find(maxPtIndices.begin(), maxPtIndices.end(), idx) == maxPtIndices.end() ){
119  (*jets)[idx].p4()+(*leps)[ 0 ].p4();
120  if( maxPt<0. || maxPt<sum.pt() ){
121  maxPt=sum.pt();
122  lepB=idx;
123  }
124  }
125  }
126 
127  // -----------------------------------------------------
128  // add jets
129  // -----------------------------------------------------
130  if( isValid(closestToWMassIndices[0], jets) ){
131  setCandidate(jets, closestToWMassIndices[0], lightQ_, jetCorrectionLevel("wQuarkMix"));
132  match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0];
133  }
134 
135  if( isValid(closestToWMassIndices[1], jets) ){
136  setCandidate(jets, closestToWMassIndices[1], lightQBar_, jetCorrectionLevel("wQuarkMix"));
137  match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
138  }
139 
140  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){
141  // if this idx is not yet contained in the list of W mass candidates...
142  if( std::find( closestToWMassIndices.begin(), closestToWMassIndices.end(), maxPtIndices[idx]) == closestToWMassIndices.end() ){
143  // ...and if it is valid
144  if( isValid(maxPtIndices[idx], jets) ){
145  setCandidate(jets, maxPtIndices[idx], hadronicB_, jetCorrectionLevel("bQuark"));
146  match[TtSemiLepEvtPartons::HadB] = maxPtIndices[idx];
147  break; // there should be no other cadidates!
148  }
149  }
150  }
151 
152  if( isValid(lepB, jets) ){
153  setCandidate(jets, lepB, leptonicB_, jetCorrectionLevel("bQuark"));
154  match[TtSemiLepEvtPartons::LepB] = lepB;
155  }
156 
157  // -----------------------------------------------------
158  // add lepton
159  // -----------------------------------------------------
160  setCandidate(leps, 0, lepton_);
161  match[TtSemiLepEvtPartons::Lepton] = 0;
162 
163  // -----------------------------------------------------
164  // add neutrino
165  // -----------------------------------------------------
166  if(neutrinoSolutionType_ == -1)
167  setCandidate(mets, 0, neutrino_);
168  else
169  setNeutrino(mets, leps, 0, neutrinoSolutionType_);
170 }
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_
TtSemiLepHypMaxSumPtWMass(const edm::ParameterSet &)
reco::ShallowClonePtrCandidate * lightQBar_
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 find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
reco::ShallowClonePtrCandidate * neutrino_
reco::ShallowClonePtrCandidate * hadronicB_
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
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 ...