CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VBF.py
Go to the documentation of this file.
1 import math
2 from ROOT import TLorentzVector
3 from PhysicsTools.HeppyCore.utils.deltar import deltaPhi, deltaR2
4 
5 class VBF( object ):
6  '''Computes and holds VBF quantities'''
7  def __init__(self, jets, diLepton, vbfMvaCalc, cjvPtCut):
8  '''jets: jets cleaned from the diLepton legs.
9  diLepton: the di-tau, for example. Necessary to compute input variables for MVA selection
10  '''
11  self.cjvPtCut = cjvPtCut
12  self.vbfMvaCalc = vbfMvaCalc
13  self.jets = jets
14  # the MET is taken from the di-lepton, because it can depend on it
15  # e.g. recoil corrections, mva met
16  self.met = diLepton.met()
17  self.leadJets = jets[:2] # taking first 2 jets (leading ones)
18  self.otherJets = jets[2:]
19  self.centralJets = self.findCentralJets( self.leadJets, self.otherJets )
20 
21  # delta eta
22  self.deta = self.leadJets[0].eta() - self.leadJets[1].eta()
23 
24  # below, the variables for the MVA selection
25  # delta phi
26  self.dphi = deltaPhi(self.leadJets[0].phi(), self.leadJets[1].phi())
27  dijetp4 = self.leadJets[0].p4() + self.leadJets[1].p4()
28  # mass of the di-jet system
29  self.mjj = dijetp4.M()
30  # pt of di-jet system
31  self.dijetpt = dijetp4.pt()
32  # phi of di-jet system
33  self.dijetphi = dijetp4.phi()
34  # higgs momentum (defined as the di-lepton momentum + the met momentum)
35  # don't access longitudinal quantities!
36  self.higgsp4 = diLepton.p4() + self.met.p4()
37  # delta phi between dijet system and higgs system
38  self.dphidijethiggs = deltaPhi( self.dijetphi, self.higgsp4.phi() )
39  # ?
40  visDiLepton = diLepton.leg1 ().p4 () + diLepton.leg2 ().p4 ()
41  self.visjeteta = min (
42  abs (self.leadJets[0].eta () - visDiLepton.eta ()),
43  abs (self.leadJets[1].eta () - visDiLepton.eta ()))
44  # visible higgs pt = di-lepton pt
45  self.ptvis = visDiLepton.pt()
46  ## self.ptvis = diLepton.pt()
47  # new VBF MVA, based on 4 variables
48  if self.vbfMvaCalc is not None:
49  self.mva = self.vbfMvaCalc.val( self.mjj,
50  abs(self.deta),
51  self.visjeteta,
52  self.ptvis )
53  else:
54  self.mva = -99.
55 
56 # double mjj , // the invariant mass of the two tag jets
57 # double dEta , // the pseudorapidity difference between the two tag jets
58 # double dPhi , // the phi difference between the two tag jets
59 # double ditau_pt , // the vector sum of the pT of the tau + electron/muon + MET
60 # double dijet_pt , // the vector sum of the pT of the two tag jets
61 # double dPhi_hj , // the phi difference between the di-tau vector and the di-jet vector
62 # double C1 , // the pseudorapidity difference between the *visible* di-tau vector and the closest tag jet
63 # double C2 // the *visible* pT of the di-tau
64 
65 
66  def findCentralJets( self, leadJets, otherJets ):
67  '''Finds all jets between the 2 leading jets, for central jet veto.'''
68  if not len(otherJets):
69  return []
70  etamin = leadJets[0].eta()
71  etamax = leadJets[1].eta()
72  if etamin > etamax:
73  etamin, etamax = etamax, etamin
74  def isCentral( jet ):
75  if jet.pt()<self.cjvPtCut:
76  return False
77  eta = jet.eta()
78  if etamin < eta and eta < etamax:
79  return True
80  else:
81  return False
82  centralJets = filter( isCentral, otherJets )
83  return centralJets
84 
85  def calcP4(self, jets):
86  '''returns the sum p4 of a collection of objects.
87  FIXME: remove this function, which is a bit stupid
88  '''
89  p4 = TLorentzVector()
90  for jet in jets:
91  p4 += TLorentzVector(jet.px(), jet.py(), jet.pz(), jet.energy())
92  return p4
93 
94  def __str__(self):
95  header = 'VBF : deta={deta:4.2f}, Mjj={mjj:4.2f}, #centjets={ncentjets}'
96  header = header.format( deta=self.deta, mjj=self.mjj, ncentjets=len(self.centralJets))
97  leadJets = map( str, self.leadJets )
98  centralJets = map( str, self.centralJets)
99  tmp = [header]
100  tmp.append('MVA input variables: dphi={dphi:4.2f}, dijetpt={dijetpt:4.2f}, dijetphi={dijetphi:4.2f}, dphidijethiggs={dphidijethiggs:4.2f}, visjeteta={visjeteta:4.2f}, ptvis={ptvis:4.2f}'.format(
101  dphi = self.dphi,
102  dijetpt = self.dijetpt,
103  dijetphi = self.dijetphi,
104  dphidijethiggs = self.dphidijethiggs,
105  visjeteta = self.visjeteta,
106  ptvis = self.ptvis
107  ))
108  tmp.append('Leading Jets:')
109  tmp.extend( leadJets )
110  tmp.append('Central Jets:')
111  tmp.extend( centralJets )
112  return '\n'.join( tmp )
dijetphi
Definition: VBF.py:33
jets
Definition: VBF.py:13
def __init__
Definition: VBF.py:7
higgsp4
Definition: VBF.py:36
dijetpt
Definition: VBF.py:31
dphi
Definition: VBF.py:26
cjvPtCut
Definition: VBF.py:11
def calcP4
Definition: VBF.py:85
def findCentralJets
Definition: VBF.py:66
def __str__
Definition: VBF.py:94
deta
Definition: VBF.py:22
double p4[4]
Definition: TauolaWrapper.h:92
otherJets
Definition: VBF.py:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
vbfMvaCalc
Definition: VBF.py:12
ptvis
Definition: VBF.py:45
leadJets
Definition: VBF.py:17
mjj
Definition: VBF.py:29
centralJets
Definition: VBF.py:19
Definition: VBF.py:5
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
mva
self.ptvis = diLepton.pt() new VBF MVA, based on 4 variables
Definition: VBF.py:49
list object
Definition: dbtoconf.py:77
dphidijethiggs
Definition: VBF.py:38
met
Definition: VBF.py:16
visjeteta
Definition: VBF.py:41