CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes
TtDilepLRSignalSelObservables Class Reference

#include <TtDilepLRSignalSelObservables.h>

Public Types

typedef std::pair< unsigned
int, bool > 
IntBoolPair
 

Public Member Functions

void jetSource (const edm::InputTag &jetSource)
 
std::vector< IntBoolPairoperator() (TtDilepEvtSolution &, const edm::Event &iEvent, bool matchOnly=false)
 
 TtDilepLRSignalSelObservables ()
 
 ~TtDilepLRSignalSelObservables ()
 

Private Types

typedef std::pair< unsigned
int, double > 
IntDblPair
 

Private Member Functions

double delta (double phi1, double phi2)
 
void fillMinMax (double v1, double v2, int obsNbr, std::vector< IntDblPair > &varList, bool match1, bool match2, std::vector< IntBoolPair > &matchList)
 

Private Attributes

int count1
 
int count2
 
int count3
 
int count4
 
int count5
 
std::vector< IntBoolPairevtselectVarMatch
 
std::vector< IntDblPairevtselectVarVal
 
edm::InputTag jetSource_
 

Detailed Description

Definition at line 15 of file TtDilepLRSignalSelObservables.h.

Member Typedef Documentation

typedef std::pair<unsigned int,bool> TtDilepLRSignalSelObservables::IntBoolPair

Definition at line 22 of file TtDilepLRSignalSelObservables.h.

typedef std::pair<unsigned int,double> TtDilepLRSignalSelObservables::IntDblPair
private

Definition at line 29 of file TtDilepLRSignalSelObservables.h.

Constructor & Destructor Documentation

TtDilepLRSignalSelObservables::TtDilepLRSignalSelObservables ( )
TtDilepLRSignalSelObservables::~TtDilepLRSignalSelObservables ( )

Definition at line 17 of file TtDilepLRSignalSelObservables.cc.

17  {
18 }

Member Function Documentation

double TtDilepLRSignalSelObservables::delta ( double  phi1,
double  phi2 
)
private

Definition at line 208 of file TtDilepLRSignalSelObservables.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi, and M_PI.

Referenced by operator()().

209 {
210  double deltaPhi = phi1 - phi2;
211  while (deltaPhi > M_PI) deltaPhi -= 2*M_PI;
212  while (deltaPhi <= -M_PI) deltaPhi += 2*M_PI;
213  return deltaPhi;
214 }
#define M_PI
Definition: BFit3D.cc:3
void TtDilepLRSignalSelObservables::fillMinMax ( double  v1,
double  v2,
int  obsNbr,
std::vector< IntDblPair > &  varList,
bool  match1,
bool  match2,
std::vector< IntBoolPair > &  matchList 
)
private

Definition at line 191 of file TtDilepLRSignalSelObservables.cc.

Referenced by operator()().

193 {
194  if (v1<v2) {
195  varList.push_back(IntDblPair(obsNbr, v1));
196  varList.push_back(IntDblPair(obsNbr+1, v2));
197  matchList.push_back(IntBoolPair(obsNbr, match1));
198  matchList.push_back(IntBoolPair(obsNbr+1, match2));
199 
200  } else {
201  varList.push_back(IntDblPair(obsNbr, v2));
202  varList.push_back(IntDblPair(obsNbr+1, v1));
203  matchList.push_back(IntBoolPair(obsNbr, match2));
204  matchList.push_back(IntBoolPair(obsNbr+1, match1));
205  }
206 }
std::pair< unsigned int, bool > IntBoolPair
std::pair< unsigned int, double > IntDblPair
void TtDilepLRSignalSelObservables::jetSource ( const edm::InputTag jetSource)
inline

Definition at line 25 of file TtDilepLRSignalSelObservables.h.

References jetSource(), and jetSource_.

Referenced by jetSource().

std::vector< TtDilepLRSignalSelObservables::IntBoolPair > TtDilepLRSignalSelObservables::operator() ( TtDilepEvtSolution solution,
const edm::Event iEvent,
bool  matchOnly = false 
)

Definition at line 21 of file TtDilepLRSignalSelObservables.cc.

References abs, count1, count2, count3, count4, count5, gather_cfg::cout, delta(), SiPixelRawToDigiRegional_cfi::deltaPhi, reco::LeafCandidate::et(), evtselectVarMatch, evtselectVarVal, fillMinMax(), MCTruth::genEvent, edm::Event::getByLabel(), TtDilepEvtSolution::getCalJetB(), TtDilepEvtSolution::getCalJetBbar(), TtDilepEvtSolution::getGenLepm(), TtDilepEvtSolution::getGenLepp(), TtDilepEvtSolution::getJetB(), TtDilepEvtSolution::getJetBbar(), TtDilepEvtSolution::getLeptNeg(), TtDilepEvtSolution::getLeptPos(), TtDilepEvtSolution::getWmDecay(), TtDilepEvtSolution::getWpDecay(), i, fwrapper::jets, jetSource_, M_PI, scaleCards::mass, reco::Particle::p4(), reco::LeafCandidate::p4(), pat::Jet::partonFlavour(), reco::LeafCandidate::pdgId(), createTree::pp, and TtDilepEvtSolution::setLRSignalEvtObservables().

23 {
24  evtselectVarVal.clear();
25  evtselectVarMatch.clear();
26 
27  // Check whether the objects are matched:
28  bool matchB1 = false;
29  bool matchB2 = false;
30  bool matchB = false;
31  bool matchBbar = false;
32  bool matchLeptPos = false;
33  bool matchLeptNeg = false;
34 
35  try {
36  // std::cout <<std::endl;
37  double dr, dr1, dr2;
38 
40  iEvent.getByLabel ("genEvt",genEvent);
41 
42  if (genEvent->isFullLeptonic()) {
43  // Match the leptons, by type and deltaR
44  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptPos(), *(solution.getGenLepp()));
45  matchLeptPos = (
46  ( ((solution.getWpDecay()=="electron")&&(std::abs(solution.getGenLepp()->pdgId())==11))
47  || ((solution.getWpDecay()=="muon")&&(std::abs(solution.getGenLepp()->pdgId())==13)) )
48  && (dr < 0.1) );
49 
50  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptNeg(), *(solution.getGenLepm()));
51  matchLeptNeg = (
52  ( ((solution.getWmDecay()=="electron")&&(std::abs(solution.getGenLepm()->pdgId())==11))
53  || ((solution.getWmDecay()=="muon")&&(std::abs(solution.getGenLepm()->pdgId())==13)) )
54  && (dr < 0.1) );
55  }
56 
57  if (genEvent->isSemiLeptonic()) {
58  int id = genEvent->singleLepton()->pdgId();
59 
60  if (id>0) {
61  dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptNeg(), *(genEvent->singleLepton()));
62  matchLeptNeg = (
63  ( ((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 = (
69  ( ((solution.getWpDecay()=="electron")&& (id==-11))
70  || ((solution.getWpDecay()=="muon") && (id==-13)) )
71  && (dr < 0.1) );
72  }
73  }
74 
75  if (genEvent->isTtBar() && genEvent->numberOfBQuarks()>1) {
76  if (solution.getJetB().partonFlavour()==5) ++count1;
77  if (solution.getJetBbar().partonFlavour()==5) ++count1;
78 
79  dr1 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetB(), *(genEvent->b()));
80  dr2 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetB(), *(genEvent->bBar()));
81 
82  matchB1= ( (dr1<0.4) || (dr2<0.4));
83  matchB = ( (solution.getJetB().partonFlavour()==5) && (dr1<0.4) );
84  if (matchB) ++count3;
85  matchB = ( (dr1<0.4) );
86  if (dr1<0.5) ++count2;
87  if (dr1<0.4) ++count4;
88  if (dr1<0.3) ++count5;
89 
90  dr1 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetBbar(), *(genEvent->b()));
91  dr2 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetBbar(), *(genEvent->bBar()));
92 
93  matchBbar = ( (solution.getJetBbar().partonFlavour()==5) && (dr2<0.4) );
94  if (matchBbar) ++count3;
95  matchBbar = ( (dr2<0.4) );
96  matchB2 = ( (dr1<0.4) || (dr2<0.4));
97  if (dr2<0.5) ++count2;
98  if (dr2<0.4) ++count4;
99  if (dr2<0.3) ++count5;
100  }
101 
102  } catch (...){std::cout << "Exception\n";}
103 
105  iEvent.getByLabel(jetSource_, jets);
106 
107  // Lower / Higher of both jet angles
108 
109  double v1 = std::abs( solution.getJetB().p4().theta() - M_PI/2 );
110  double v2 = std::abs( solution.getJetBbar().p4().theta() - M_PI/2 ) ;
111  fillMinMax(v1, v2, 1, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
112 
113  // Lower / Higher of both jet pT
114 
115  double pt1 = solution.getJetB().p4().pt();
116  double pt2 = solution.getJetBbar().p4().pt();
117  fillMinMax(pt1, pt2, 3, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
118 
119  // Lower / Higher of both lepton pT
120 
121  pt1 = solution.getLeptPos().p4().pt();
122  pt2 = solution.getLeptNeg().p4().pt();
123  fillMinMax(pt1, pt2, 5, evtselectVarVal, matchLeptPos, matchLeptNeg, evtselectVarMatch);
124 
125  // delta theta btw the b-jets
126 
127  double deltaPhi = std::abs( delta(solution.getJetB().p4().phi(),
128  solution.getJetBbar().p4().phi()) );
129 
130  evtselectVarVal.push_back(IntDblPair(7, deltaPhi));
131  evtselectVarMatch.push_back(IntBoolPair(7, matchB1&&matchB2));
132 
133  // delta phi btw the b-jets
134 
135  double deltaTheta = std::abs( delta (solution.getJetBbar().p4().theta(),
136  solution.getJetB().p4().theta() ) );
137 
138  evtselectVarVal.push_back(IntDblPair(8, deltaTheta));
139  evtselectVarMatch.push_back(IntBoolPair(8, matchB1&&matchB2));
140 
141  // Lower / Higher of phi difference between the b and associated lepton
142 
143  double deltaPhi1 = std::abs ( delta( solution.getJetB().p4().phi(),
144  solution.getLeptPos().p4().phi() ) );
145  double deltaPhi2 = std::abs ( delta( solution.getJetBbar().p4().phi(),
146  solution.getLeptNeg().p4().phi() ) );
147 
148  fillMinMax(deltaPhi1, deltaPhi2, 9, evtselectVarVal,
149  matchB&&matchLeptPos, matchBbar&&matchLeptNeg, evtselectVarMatch);
150 
151  // Lower / Higher of theta difference between the b and associated lepton
152 
153  double deltaTheta1 = std::abs( solution.getJetB().p4().theta() -
154  solution.getLeptPos().p4().theta() );
155  double deltaTheta2 = std::abs( solution.getJetBbar().p4().theta() -
156  solution.getLeptNeg().p4().theta() );
157  fillMinMax(deltaTheta1, deltaTheta2, 11, evtselectVarVal,
158  matchB&&matchLeptPos, matchBbar&&matchLeptNeg, evtselectVarMatch);
159 
160  // Invariant Mass of lepton pair
161 
162  math::XYZTLorentzVector pp = solution.getLeptPos().p4() + solution.getLeptNeg().p4();
163  double mass = pp.mass();
164  evtselectVarVal.push_back(IntDblPair(13, mass));
165  evtselectVarMatch.push_back(IntBoolPair(13, matchLeptNeg&&matchLeptPos));
166 
167  evtselectVarVal.push_back(IntDblPair(13, mass));
168  evtselectVarMatch.push_back(IntBoolPair(13, matchLeptNeg&&matchLeptPos));
169 
170  std::vector <pat::Jet> jet3;
171  for (int i=0;i<(int)jets->size();++i) {
172 if ( ((*jets)[i].et()<solution.getJetB().et()) && ((*jets)[i].et()<solution.getJetBbar().et())) {jet3.push_back((*jets)[i]);
173 }}
174  double jet1Ratio = 0., jet2Ratio = 0.;
175  if (jet3.size()>0) {
176  jet1Ratio = jet3[0].et()/solution.getJetB().et();
177  jet2Ratio = jet3[0].et()/solution.getJetBbar().et();
178  }
179  fillMinMax(jet1Ratio, jet2Ratio, 14, evtselectVarVal,
180  matchB1, matchB2, evtselectVarMatch);
181 
182  evtselectVarVal.push_back(IntDblPair(16, jets->size()));
183  evtselectVarMatch.push_back(IntBoolPair(16, matchB&&matchBbar));
184 
185 
186  if (!matchOnly) solution.setLRSignalEvtObservables(evtselectVarVal);
187  return evtselectVarMatch;
188 }
pat::Jet getCalJetBbar() const
const reco::GenParticle * getGenLepm() const
double delta(double phi1, double phi2)
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
void setLRSignalEvtObservables(std::vector< std::pair< unsigned int, double > >)
tuple pp
Definition: createTree.py:15
virtual double et() const
transverse energy
void fillMinMax(double v1, double v2, int obsNbr, std::vector< IntDblPair > &varList, bool match1, bool match2, std::vector< IntBoolPair > &matchList)
Definition: deltaR.h:51
std::pair< unsigned int, bool > IntBoolPair
#define abs(x)
Definition: mlp_lapack.h:159
const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: Particle.h:64
std::string getWmDecay() const
std::string getWpDecay() const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
std::pair< unsigned int, double > IntDblPair
reco::Particle getLeptNeg() const
vector< PseudoJet > jets
pat::Jet getJetBbar() const
tuple genEvent
Definition: MCTruth.py:33
reco::Particle getLeptPos() const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
int partonFlavour() const
return the parton-based flavour of the jet
Definition: Jet.cc:196
std::vector< IntBoolPair > evtselectVarMatch
#define M_PI
Definition: BFit3D.cc:3
pat::Jet getJetB() const
pat::Jet getCalJetB() const
const reco::GenParticle * getGenLepp() const
tuple mass
Definition: scaleCards.py:27
tuple cout
Definition: gather_cfg.py:121
virtual const LorentzVector & p4() const
four-momentum Lorentz vector

Member Data Documentation

int TtDilepLRSignalSelObservables::count1
private

Definition at line 40 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()(), and TtDilepLRSignalSelObservables().

int TtDilepLRSignalSelObservables::count2
private

Definition at line 40 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()(), and TtDilepLRSignalSelObservables().

int TtDilepLRSignalSelObservables::count3
private

Definition at line 40 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()(), and TtDilepLRSignalSelObservables().

int TtDilepLRSignalSelObservables::count4
private

Definition at line 40 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()(), and TtDilepLRSignalSelObservables().

int TtDilepLRSignalSelObservables::count5
private

Definition at line 40 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()(), and TtDilepLRSignalSelObservables().

std::vector< IntBoolPair > TtDilepLRSignalSelObservables::evtselectVarMatch
private

Definition at line 39 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()().

std::vector< IntDblPair > TtDilepLRSignalSelObservables::evtselectVarVal
private

Definition at line 38 of file TtDilepLRSignalSelObservables.h.

Referenced by operator()().

edm::InputTag TtDilepLRSignalSelObservables::jetSource_
private

Definition at line 36 of file TtDilepLRSignalSelObservables.h.

Referenced by jetSource(), and operator()().