CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepJetComb.cc
Go to the documentation of this file.
3 
4 #include "Math/VectorUtil.h"
5 
7 {
8 }
9 
10 TtSemiLepJetComb::TtSemiLepJetComb(const std::vector<pat::Jet>& jets, const std::vector<int>& combination,
11  const math::XYZTLorentzVector& lepton, const pat::MET& neutrino)
12 {
13  // receive right jet association
14  // from jet-parton matching
15  hadQJet_ = jets[ combination[TtSemiLepEvtPartons::LightQ ] ];
16  hadQBarJet_ = jets[ combination[TtSemiLepEvtPartons::LightQBar] ];
17  hadBJet_ = jets[ combination[TtSemiLepEvtPartons::HadB ] ];
18  lepBJet_ = jets[ combination[TtSemiLepEvtPartons::LepB ] ];
19  lepton_ = lepton;
20  neutrino_ = neutrino;
21  // create mother candidates from
22  // final-state candidates
23  deduceMothers();
24 }
25 
27 {
28 }
29 
31 {
32  switch(var){
33  case JetComb::kMass : return top(decay).mass();
34  case JetComb::kPt : return top(decay).pt();
35  case JetComb::kEta : return top(decay).eta();
36  case JetComb::kPhi : return top(decay).phi();
37  case JetComb::kTheta : return top(decay).theta();
38  };
39  return -9999.;
40 }
41 
43 {
44  switch(var){
45  case JetComb::kMass : return wBoson(decay).mass();
46  case JetComb::kPt : return wBoson(decay).pt();
47  case JetComb::kEta : return wBoson(decay).eta();
48  case JetComb::kPhi : return wBoson(decay).phi();
49  case JetComb::kTheta : return wBoson(decay).theta();
50  };
51  return -9999.;
52 }
53 
55 {
56  switch(var){
57  case JetComb::kMass : return -9999.;
58  case JetComb::kPt : return bQuark(decay).p4().pt();
59  case JetComb::kEta : return bQuark(decay).p4().eta();
60  case JetComb::kPhi : return bQuark(decay).p4().phi();
61  case JetComb::kTheta : return bQuark(decay).p4().theta();
62  };
63  return -9999.;
64 }
65 
67 {
68  switch(var){
69  case JetComb::kMass : return -9999.;
70  case JetComb::kPt : return lightQ(qbar).p4().pt();
71  case JetComb::kEta : return lightQ(qbar).p4().eta();
72  case JetComb::kPhi : return lightQ(qbar).p4().phi();
73  case JetComb::kTheta : return lightQ(qbar).p4().theta();
74  };
75  return -9999.;
76 }
77 
79 {
80  switch(var){
81  case JetComb::kMass : return -9999.;
82  case JetComb::kPt : return lepton_.pt();
83  case JetComb::kEta : return lepton_.eta();
84  case JetComb::kPhi : return lepton_.phi();
85  case JetComb::kTheta : return lepton_.theta();
86  };
87  return -9999.;
88 }
89 
91 {
92  switch(var){
93  case JetComb::kMass : return -9999.;
94  case JetComb::kPt : return neutrino_.p4().pt();
95  case JetComb::kEta : return neutrino_.p4().eta();
96  case JetComb::kPhi : return neutrino_.p4().phi();
97  case JetComb::kTheta : return neutrino_.p4().theta();
98  };
99  return -9999.;
100 }
101 
103 {
104  switch(comp){
105  case JetComb::kDeltaM : return top(JetComb::kHad).mass() - top(JetComb::kLep).mass();
108  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (top(JetComb::kHad), top(JetComb::kLep));
109  };
110  return -9999.;
111 }
112 
114 {
115  switch(comp){
116  case JetComb::kDeltaM : return wBoson(JetComb::kHad).mass() - wBoson(JetComb::kLep).mass();
119  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (wBoson(JetComb::kHad), wBoson(JetComb::kLep));
120  };
121  return -9999.;
122 }
123 
125 {
126  switch(comp){
127  case JetComb::kDeltaM : return -9999.;
130  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (bQuark(JetComb::kHad).p4(), bQuark(JetComb::kLep).p4());
131  };
132  return -9999.;
133 }
134 
136 {
137  switch(comp){
138  case JetComb::kDeltaM : return -9999.;
139  case JetComb::kDeltaR : return ROOT::Math::VectorUtil::DeltaR (lightQ().p4(), lightQ(true).p4());
140  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(lightQ().p4(), lightQ(true).p4()));
141  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (lightQ().p4(), lightQ(true).p4());
142  };
143  return -9999.;
144 }
145 
147 {
148  switch(comp){
149  case JetComb::kDeltaM : return -9999.;
152  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (lepton_, neutrino_.p4());
153  };
154  return -9999.;
155 }
156 
158 {
159  switch(comp){
160  case JetComb::kDeltaM : return top(dec1).mass() - wBoson(dec2).mass();
161  case JetComb::kDeltaR : return ROOT::Math::VectorUtil::DeltaR (top(dec1), wBoson(dec2));
162  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(top(dec1), wBoson(dec2)));
163  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (top(dec1), wBoson(dec2));
164  };
165  return -9999.;
166 }
167 
169 {
170  switch(comp){
171  case JetComb::kDeltaM : return -9999.;
172  case JetComb::kDeltaR : return ROOT::Math::VectorUtil::DeltaR (top(dec1), bQuark(dec2).p4());
173  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(top(dec1), bQuark(dec2).p4()));
174  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (top(dec1), bQuark(dec2).p4());
175  };
176  return -9999.;
177 }
178 
180 {
181  switch(comp){
182  case JetComb::kDeltaM : return -9999.;
183  case JetComb::kDeltaR : return ROOT::Math::VectorUtil::DeltaR (wBoson(dec1), bQuark(dec2).p4());
184  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(wBoson(dec1), bQuark(dec2).p4()));
185  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (wBoson(dec1), bQuark(dec2).p4());
186  };
187  return -9999.;
188 }
189 
191 {
192  switch(comp){
193  case JetComb::kDeltaM : return -9999.;
195  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(top(decay), lepton_));
196  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (top(decay), lepton_);
197  };
198  return -9999.;
199 }
200 
202 {
203  switch(comp){
204  case JetComb::kDeltaM : return -9999.;
206  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(top(decay), neutrino_.p4()));
207  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (top(decay), neutrino_.p4());
208  };
209  return -9999.;
210 }
211 
213 {
214  switch(comp){
215  case JetComb::kDeltaM : return -9999.;
218  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (wBoson(decay), lepton_);
219  };
220  return -9999.;
221 }
222 
224 {
225  switch(comp){
226  case JetComb::kDeltaM : return -9999.;
229  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (wBoson(decay), neutrino_.p4());
230  };
231  return -9999.;
232 }
233 
235 {
236  switch(comp){
237  case JetComb::kDeltaM : return -9999.;
240  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (bQuark(decay).p4(), lepton_);
241  };
242  return -9999.;
243 }
244 
246 {
247  switch(comp){
248  case JetComb::kDeltaM : return -9999.;
250  case JetComb::kDeltaPhi : return fabs(ROOT::Math::VectorUtil::DeltaPhi(bQuark(decay).p4(), neutrino_.p4()));
251  case JetComb::kDeltaTheta : return ROOT::Math::VectorUtil::Angle (bQuark(decay).p4(), neutrino_.p4());
252  };
253  return -9999.;
254 }
255 
257 {
258  return top(JetComb::kHad).pt()/(top(JetComb::kHad).pt() +
259  (lightQ() .p4()+lightQ(true).p4()+bQuark(JetComb::kLep).p4()).pt() +
262  );
263 }
264 
266 {
267  return (bQuark(JetComb::kHad).p4().pt()+bQuark(JetComb::kLep).p4().pt())/(lightQ().p4().pt()+lightQ(true).p4().pt());
268 }
269 
271 {
272  switch(op){
273  case JetComb::kAdd : return bTag(JetComb::kHad, algo) + bTag(JetComb::kLep, algo);
274  case JetComb::kMult : return bTag(JetComb::kHad, algo) * bTag(JetComb::kLep, algo);
275  };
276  return -9999.;
277 }
278 
280 {
281  switch(op){
282  case JetComb::kAdd : return bTag(lightQ(), algo) + bTag(lightQ(true), algo);
283  case JetComb::kMult : return bTag(lightQ(), algo) * bTag(lightQ(true), algo);
284  };
285  return -9999.;
286 }
287 
288 // ----------------------------------------------------------------------
289 // private methods
290 // ----------------------------------------------------------------------
291 
292 void
294 {
295  hadW_ = hadQJet_.p4() + hadQBarJet_.p4();
296  lepW_ = lepton_ + neutrino_ .p4();
297  hadTop_ = hadW_ + hadBJet_ .p4();
298  lepTop_ = lepW_ + lepBJet_ .p4();
299 }
300 
302 {
303  switch(algo){
304  case JetComb::kTrackCountHighEff : return jet.bDiscriminator("trackCountingHighEffBJetTags" );
305  case JetComb::kTrackCountHighPur : return jet.bDiscriminator("trackCountingHighPurBJetTags" );
306  case JetComb::kSoftMuon : return jet.bDiscriminator("softMuonBJetTags" );
307  case JetComb::kSoftMuonByPt : return jet.bDiscriminator("softMuonByPtBJetTags" );
308  case JetComb::kSofMuonByIP3d : return jet.bDiscriminator("softMuonByIP3dBJetTags" );
309  case JetComb::kSoftElec : return jet.bDiscriminator("softElectronBJetTags" );
310  case JetComb::kProbability : return jet.bDiscriminator("jetProbabilityBJetTags" );
311  case JetComb::kBProbability : return jet.bDiscriminator("jetBProbabilityBJetTags" );
312  case JetComb::kSimpleSecondVtx : return jet.bDiscriminator("simpleSecondaryVertexBJetTags" );
313  case JetComb::kCombSecondVtx : return jet.bDiscriminator("combinedSecondaryVertexBJetTags" );
314  case JetComb::kCombSecondVtxMVA : return jet.bDiscriminator("combinedSecondaryVertexMVABJetTags");
315  };
316  return -9999.;
317 }
Analysis-level MET class.
Definition: MET.h:43
double relativePtHadronicTop() const
math::XYZTLorentzVector lepton_
lepton 4-vector
math::XYZTLorentzVector hadTop_
hadronic top 4-vector
CompType
supported comparison types
math::XYZTLorentzVector lepW_
leptonic W 4-vector
double combinedBTags(JetComb::BTagAlgo algo, JetComb::Operator op) const
combined b-tag values of the two b candidates
math::XYZTLorentzVector lepTop_
leptonic top 4-vector
double compareWB(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const
comparison between the W and the b candidate
double neutrinoVar(JetComb::VarType var) const
neutrino candidate variable
double compareHadBLepB(JetComb::CompType comp) const
comparison between the two b candidates
~TtSemiLepJetComb()
default destructor
double compareHadWLepW(JetComb::CompType comp) const
comparison between the two W candidates
TtSemiLepJetComb()
emtpy constructor
double bQuarkVar(JetComb::DecayType decay, JetComb::VarType var) const
b quark candidate variable
float bDiscriminator(const std::string &theLabel) const
-— methods for accessing b-tagging info -—
const math::XYZTLorentzVector & top(JetComb::DecayType decay) const
return leptonic or hadronic top candidate depending on argument
double compareHadTopLepTop(JetComb::CompType comp) const
comparison between the two top candidates
double compareWLepton(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the W and the lepton candidate
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double compareBLepton(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the b and the lepton candidate
double lightQVar(bool qbar, JetComb::VarType var) const
light quark candidate variable
double compareBNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the b and the neutrino candidate
double compareLightQuarks(JetComb::CompType comp) const
comparison between the lightQ and the lightQBar candidate
double compareLeptonNeutrino(JetComb::CompType comp) const
comparison between the lepton and the neutrino candidate
double p4[4]
Definition: TauolaWrapper.h:92
vector< PseudoJet > jets
Operator
operators for combining variables
double compareTopLepton(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the top and the lepton candidate
double compareTopW(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const
comparison between the top and the W candidate
double compareTopB(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const
comparison between the top and the b candidate
const pat::Jet & lightQ(bool qbar=false) const
return lightQ or lightQBar candidate depending on argument
pat::Jet hadBJet_
hadronic b jet
double topVar(JetComb::DecayType decay, JetComb::VarType var) const
top quark candidate variable
double leptonVar(JetComb::VarType var) const
lepton candidate variable
double wBosonVar(JetComb::DecayType decay, JetComb::VarType var) const
W boson candidate variable.
double combinedBTagsForLightQuarks(JetComb::BTagAlgo algo, JetComb::Operator op) const
combined b-tag values of the two light quark candidates
double compareWNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the W and the neutrino candidate
pat::Jet hadQJet_
lightQ jet
pat::Jet lepBJet_
leptonic b jet
math::XYZTLorentzVector hadW_
hadronic W 4-vector
Analysis-level calorimeter jet class.
Definition: Jet.h:77
VarType
supported std single variable types
pat::Jet hadQBarJet_
lightQBar jet
BTagAlgo
b-tagging algorithms
const math::XYZTLorentzVector & wBoson(JetComb::DecayType decay) const
return leptonic or hadronic W candidate depending on argument
double compareTopNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the top and the neutrino candidate
void deduceMothers()
reconstruct mother candidates from final state candidates
const pat::Jet & bQuark(JetComb::DecayType decay) const
return leptonic or hadronic b candidate depending on argument
double bTag(JetComb::DecayType decay, JetComb::BTagAlgo algo) const
b-tag value of a b candidate
double bOverLightQPt() const
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
pat::MET neutrino_
neutrino candidate