CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtDilepLRSignalSelObservables.cc
Go to the documentation of this file.
7 
8 /************** Definition of the functions of the class ***************/
9 
10 //Constructor
12  edm::ConsumesCollector &&iC, const edm::EDGetTokenT<std::vector<pat::Jet> > &jetSourceToken)
13  : jetSourceToken_(jetSourceToken), genEvtToken_(iC.consumes<TtGenEvent>(edm::InputTag("genEvt"))) {
14  count1 = 0;
15  count2 = 0;
16  count3 = 0;
17  count4 = 0;
18  count5 = 0;
19  count3 = 0;
20 }
21 
22 // Destructor
24 
25 std::vector<TtDilepLRSignalSelObservables::IntBoolPair> TtDilepLRSignalSelObservables::operator()(
26  TtDilepEvtSolution &solution, const edm::Event &iEvent, bool matchOnly) {
27  evtselectVarVal.clear();
28  evtselectVarMatch.clear();
29 
30  // Check whether the objects are matched:
31  bool matchB1 = false;
32  bool matchB2 = false;
33  bool matchB = false;
34  bool matchBbar = false;
35  bool matchLeptPos = false;
36  bool matchLeptNeg = false;
37 
39  iEvent.getByToken(genEvtToken_, genEvent);
40 
41  if (!genEvent.failedToGet() && genEvent.isValid()) {
42  // std::cout <<std::endl;
43  double dr, dr1, dr2;
44 
45  if (genEvent->isFullLeptonic()) {
46  // Match the leptons, by type and deltaR
47  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptPos(), *(solution.getGenLepp()));
48  matchLeptPos = ((((solution.getWpDecay() == "electron") && (std::abs(solution.getGenLepp()->pdgId()) == 11)) ||
49  ((solution.getWpDecay() == "muon") && (std::abs(solution.getGenLepp()->pdgId()) == 13))) &&
50  (dr < 0.1));
51 
52  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptNeg(), *(solution.getGenLepm()));
53  matchLeptNeg = ((((solution.getWmDecay() == "electron") && (std::abs(solution.getGenLepm()->pdgId()) == 11)) ||
54  ((solution.getWmDecay() == "muon") && (std::abs(solution.getGenLepm()->pdgId()) == 13))) &&
55  (dr < 0.1));
56  }
57 
58  if (genEvent->isSemiLeptonic()) {
59  int id = genEvent->singleLepton()->pdgId();
60 
61  if (id > 0) {
62  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptNeg(), *(genEvent->singleLepton()));
63  matchLeptNeg = ((((solution.getWmDecay() == "electron") && (id == 11)) ||
64  ((solution.getWmDecay() == "muon") && (id == 13))) &&
65  (dr < 0.1));
66  } else {
67  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptPos(), *(genEvent->singleLepton()));
68  matchLeptPos = ((((solution.getWpDecay() == "electron") && (id == -11)) ||
69  ((solution.getWpDecay() == "muon") && (id == -13))) &&
70  (dr < 0.1));
71  }
72  }
73 
74  if (genEvent->isTtBar() && genEvent->numberOfBQuarks() > 1) {
75  if (solution.getJetB().partonFlavour() == 5)
76  ++count1;
77  if (solution.getJetBbar().partonFlavour() == 5)
78  ++count1;
79 
80  dr1 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetB(), *(genEvent->b()));
81  dr2 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetB(), *(genEvent->bBar()));
82 
83  matchB1 = ((dr1 < 0.4) || (dr2 < 0.4));
84  matchB = ((solution.getJetB().partonFlavour() == 5) && (dr1 < 0.4));
85  if (matchB)
86  ++count3;
87  matchB = ((dr1 < 0.4));
88  if (dr1 < 0.5)
89  ++count2;
90  if (dr1 < 0.4)
91  ++count4;
92  if (dr1 < 0.3)
93  ++count5;
94 
95  dr1 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetBbar(), *(genEvent->b()));
96  dr2 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetBbar(), *(genEvent->bBar()));
97 
98  matchBbar = ((solution.getJetBbar().partonFlavour() == 5) && (dr2 < 0.4));
99  if (matchBbar)
100  ++count3;
101  matchBbar = ((dr2 < 0.4));
102  matchB2 = ((dr1 < 0.4) || (dr2 < 0.4));
103  if (dr2 < 0.5)
104  ++count2;
105  if (dr2 < 0.4)
106  ++count4;
107  if (dr2 < 0.3)
108  ++count5;
109  }
110  }
111 
113  iEvent.getByToken(jetSourceToken_, jets);
114 
115  // Lower / Higher of both jet angles
116 
117  double v1 = std::abs(solution.getJetB().p4().theta() - M_PI / 2);
118  double v2 = std::abs(solution.getJetBbar().p4().theta() - M_PI / 2);
119  fillMinMax(v1, v2, 1, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
120 
121  // Lower / Higher of both jet pT
122 
123  double pt1 = solution.getJetB().p4().pt();
124  double pt2 = solution.getJetBbar().p4().pt();
125  fillMinMax(pt1, pt2, 3, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
126 
127  // Lower / Higher of both lepton pT
128 
129  pt1 = solution.getLeptPos().p4().pt();
130  pt2 = solution.getLeptNeg().p4().pt();
131  fillMinMax(pt1, pt2, 5, evtselectVarVal, matchLeptPos, matchLeptNeg, evtselectVarMatch);
132 
133  // delta theta btw the b-jets
134 
135  double deltaPhi = std::abs(delta(solution.getJetB().p4().phi(), solution.getJetBbar().p4().phi()));
136 
137  evtselectVarVal.push_back(IntDblPair(7, deltaPhi));
138  evtselectVarMatch.push_back(IntBoolPair(7, matchB1 && matchB2));
139 
140  // delta phi btw the b-jets
141 
142  double deltaTheta = std::abs(delta(solution.getJetBbar().p4().theta(), solution.getJetB().p4().theta()));
143 
144  evtselectVarVal.push_back(IntDblPair(8, deltaTheta));
145  evtselectVarMatch.push_back(IntBoolPair(8, matchB1 && matchB2));
146 
147  // Lower / Higher of phi difference between the b and associated lepton
148 
149  double deltaPhi1 = std::abs(delta(solution.getJetB().p4().phi(), solution.getLeptPos().p4().phi()));
150  double deltaPhi2 = std::abs(delta(solution.getJetBbar().p4().phi(), solution.getLeptNeg().p4().phi()));
151 
152  fillMinMax(
153  deltaPhi1, deltaPhi2, 9, evtselectVarVal, matchB && matchLeptPos, matchBbar && matchLeptNeg, evtselectVarMatch);
154 
155  // Lower / Higher of theta difference between the b and associated lepton
156 
157  double deltaTheta1 = std::abs(solution.getJetB().p4().theta() - solution.getLeptPos().p4().theta());
158  double deltaTheta2 = std::abs(solution.getJetBbar().p4().theta() - solution.getLeptNeg().p4().theta());
159  fillMinMax(deltaTheta1,
160  deltaTheta2,
161  11,
163  matchB && matchLeptPos,
164  matchBbar && matchLeptNeg,
166 
167  // Invariant Mass of lepton pair
168 
169  math::XYZTLorentzVector pp = solution.getLeptPos().p4() + solution.getLeptNeg().p4();
170  double mass = pp.mass();
171  evtselectVarVal.push_back(IntDblPair(13, mass));
172  evtselectVarMatch.push_back(IntBoolPair(13, matchLeptNeg && matchLeptPos));
173 
174  evtselectVarVal.push_back(IntDblPair(13, mass));
175  evtselectVarMatch.push_back(IntBoolPair(13, matchLeptNeg && matchLeptPos));
176 
177  std::vector<pat::Jet> jet3;
178  for (int i = 0; i < (int)jets->size(); ++i) {
179  if (((*jets)[i].et() < solution.getJetB().et()) && ((*jets)[i].et() < solution.getJetBbar().et())) {
180  jet3.push_back((*jets)[i]);
181  }
182  }
183  double jet1Ratio = 0., jet2Ratio = 0.;
184  if (!jet3.empty()) {
185  jet1Ratio = jet3[0].et() / solution.getJetB().et();
186  jet2Ratio = jet3[0].et() / solution.getJetBbar().et();
187  }
188  fillMinMax(jet1Ratio, jet2Ratio, 14, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
189 
190  evtselectVarVal.push_back(IntDblPair(16, jets->size()));
191  evtselectVarMatch.push_back(IntBoolPair(16, matchB && matchBbar));
192 
193  if (!matchOnly)
195  return evtselectVarMatch;
196 }
197 
199  double v2,
200  int obsNbr,
201  std::vector<IntDblPair> &varList,
202  bool match1,
203  bool match2,
204  std::vector<IntBoolPair> &matchList) {
205  if (v1 < v2) {
206  varList.push_back(IntDblPair(obsNbr, v1));
207  varList.push_back(IntDblPair(obsNbr + 1, v2));
208  matchList.push_back(IntBoolPair(obsNbr, match1));
209  matchList.push_back(IntBoolPair(obsNbr + 1, match2));
210 
211  } else {
212  varList.push_back(IntDblPair(obsNbr, v2));
213  varList.push_back(IntDblPair(obsNbr + 1, v1));
214  matchList.push_back(IntBoolPair(obsNbr, match2));
215  matchList.push_back(IntBoolPair(obsNbr + 1, match1));
216  }
217 }
218 
219 double TtDilepLRSignalSelObservables::delta(double phi1, double phi2) {
220  double deltaPhi = phi1 - phi2;
221  while (deltaPhi > M_PI)
222  deltaPhi -= 2 * M_PI;
223  while (deltaPhi <= -M_PI)
224  deltaPhi += 2 * M_PI;
225  return deltaPhi;
226 }
pat::Jet getCalJetBbar() const
const reco::GenParticle * getGenLepm() const
double delta(double phi1, double phi2)
tuple pp
Definition: createTree.py:17
void fillMinMax(double v1, double v2, int obsNbr, std::vector< IntDblPair > &varList, bool match1, bool match2, std::vector< IntBoolPair > &matchList)
Definition: deltaR.h:58
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::vector< IntBoolPair > operator()(TtDilepEvtSolution &, const edm::Event &iEvent, bool matchOnly=false)
const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: Particle.h:88
const LorentzVector & p4() const final
four-momentum Lorentz vector
std::string getWmDecay() const
std::string getWpDecay() const
int pdgId() const final
PDG identifier.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
std::pair< unsigned int, bool > IntBoolPair
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
void setLRSignalEvtObservables(const std::vector< std::pair< unsigned int, double > > &)
std::pair< unsigned int, double > IntDblPair
reco::Particle getLeptNeg() const
vector< PseudoJet > jets
pat::Jet getJetBbar() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< TtGenEvent > genEvtToken_
tuple genEvent
Definition: MCTruth.py:34
TtDilepLRSignalSelObservables(edm::ConsumesCollector &&iC, const edm::EDGetTokenT< std::vector< pat::Jet > > &jetSourceToken)
int partonFlavour() const
return the parton-based flavour of the jet
reco::Particle getLeptPos() const
bool isValid() const
Definition: HandleBase.h:70
#define M_PI
bool failedToGet() const
Definition: HandleBase.h:72
pat::Jet getJetB() const
edm::EDGetTokenT< std::vector< pat::Jet > > jetSourceToken_
pat::Jet getCalJetB() const
const reco::GenParticle * getGenLepp() const
double et() const final
transverse energy
genEvtToken_(mayConsume< TtGenEvent >(genEvt_))
std::vector< IntBoolPair > evtselectVarMatch