CMS 3D CMS Logo

BoostedTopProducer.cc
Go to the documentation of this file.
1 #include "BoostedTopProducer.h"
3 
4 #include <string>
5 #include <sstream>
6 #include <iostream>
7 
8 using std::cout;
9 using std::endl;
10 using std::string;
11 
12 //
13 // constants, enums and typedefs
14 //
15 
16 //
17 // static data member definitions
18 //
19 
20 //
21 // constructors and destructor
22 //
24  : eleToken_(consumes<std::vector<pat::Electron> >(iConfig.getParameter<edm::InputTag>("electronLabel"))),
25  muoToken_(consumes<std::vector<pat::Muon> >(iConfig.getParameter<edm::InputTag>("muonLabel"))),
26  jetToken_(consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("jetLabel"))),
27  metToken_(consumes<std::vector<pat::MET> >(iConfig.getParameter<edm::InputTag>("metLabel"))),
28  solToken_(mayConsume<TtSemiLeptonicEvent>(iConfig.getParameter<edm::InputTag>("solLabel"))),
29  caloIsoCut_(iConfig.getParameter<double>("caloIsoCut")),
30  mTop_(iConfig.getParameter<double>("mTop")) {
31  //register products
32  produces<std::vector<reco::CompositeCandidate> >();
33 }
34 
36 
37 //
38 // member functions
39 //
40 
41 // ------------ method called to produce the data ------------
43  using namespace edm;
44 
45  bool debug = false;
46 
47  // -----------------------------------------------------
48  // get the bare PAT objects
49  // -----------------------------------------------------
51  iEvent.getByToken(muoToken_, muonHandle);
52  std::vector<pat::Muon> const& muons = *muonHandle;
53 
55  iEvent.getByToken(jetToken_, jetHandle);
56  std::vector<pat::Jet> const& jets = *jetHandle;
57 
59  iEvent.getByToken(eleToken_, electronHandle);
60  std::vector<pat::Electron> const& electrons = *electronHandle;
61 
63  iEvent.getByToken(metToken_, metHandle);
64  std::vector<pat::MET> const& mets = *metHandle;
65 
66  // -----------------------------------------------------
67  // Event Preselection:
68  // <= 1 isolated electron or muon
69  // >= 1 electron or muon
70  // >= 2 jets
71  // >= 1 missing et
72  //
73  // To explain:
74  // We want to look at leptons within "top jets" in some
75  // cases. This means the isolation will kill those events.
76  // However, if there IS an isolated lepton, we want only
77  // one of them.
78  //
79  // So to select the prompt W lepton, the logic is:
80  // 1. If there is an isolated lepton, accept it as the W lepton.
81  // 2. Else, take the highest Pt lepton (possibly non-isolated)
82  //
83  // -----------------------------------------------------
84  bool preselection = true;
85 
86  // This will hold the prompt W lepton candidate, and a
87  // maximum pt decision variable
88  double maxWLeptonPt = -1;
89  //reco::Candidate const * Wlepton = 0;
90 
91  // ----------------------
92  // Find isolated muons, and highest pt lepton
93  // ----------------------
94  std::vector<pat::Muon>::const_iterator isolatedMuon = muons.end();
95  std::vector<pat::Muon>::const_iterator muon = muons.end();
96  unsigned int nIsolatedMuons = 0;
97  std::vector<pat::Muon>::const_iterator muonIt = muons.begin(), muonEnd = muons.end();
98  for (; muonIt != muonEnd; ++muonIt) {
99  // Find highest pt lepton
100  double pt = muonIt->pt();
101  if (pt > maxWLeptonPt) {
102  maxWLeptonPt = pt;
103  muon = muonIt;
104  }
105 
106  // Find any isolated muons
107  double caloIso = muonIt->caloIso();
108  if (caloIso >= 0 && caloIso < caloIsoCut_) {
109  nIsolatedMuons++;
110  isolatedMuon = muonIt;
111  }
112  }
113 
114  // ----------------------
115  // Find isolated electrons, and highest pt lepton
116  // ----------------------
117  std::vector<pat::Electron>::const_iterator isolatedElectron = electrons.end();
118  std::vector<pat::Electron>::const_iterator electron = electrons.end();
119  unsigned int nIsolatedElectrons = 0;
120  std::vector<pat::Electron>::const_iterator electronIt = electrons.begin(), electronEnd = electrons.end();
121  for (; electronIt != electronEnd; ++electronIt) {
122  // Find highest pt lepton
123  double pt = electronIt->pt();
124  if (pt > maxWLeptonPt) {
125  maxWLeptonPt = pt;
126  electron = electronIt;
127  }
128 
129  // Find any isolated electrons
130  double caloIso = electronIt->caloIso();
131  if (caloIso >= 0 && caloIso < caloIsoCut_) {
132  nIsolatedElectrons++;
133  isolatedElectron = electronIt;
134  }
135  }
136 
137  // ----------------------
138  // Now decide on the "prompt" lepton from the W:
139  // Choose isolated leptons over all, and if no isolated,
140  // then take highest pt lepton.
141  // ----------------------
142  bool isMuon = true;
143  if (isolatedMuon != muonEnd) {
144  muon = isolatedMuon;
145  isMuon = true;
146  } else if (isolatedElectron != electronEnd) {
147  electron = isolatedElectron;
148  isMuon = false;
149  } else {
150  // Set to the highest pt lepton
151  if (muon != muonEnd && electron == electronEnd)
152  isMuon = true;
153  else if (muon == muonEnd && electron != electronEnd)
154  isMuon = false;
155  else if (muon != muonEnd && electron != electronEnd) {
156  isMuon = muon->pt() > electron->pt();
157  }
158  }
159 
160  // ----------------------
161  // Veto events that have more than one isolated lepton
162  // ----------------------
163  int nIsolatedLeptons = nIsolatedMuons + nIsolatedElectrons;
164  if (nIsolatedLeptons > 1) {
165  preselection = false;
166  }
167 
168  // ----------------------
169  // Veto events that have no prompt lepton candidates
170  // ----------------------
171  if (muon == muonEnd && electron == electronEnd) {
172  preselection = false;
173  }
174 
175  // ----------------------
176  // Veto events with < 2 jets or no missing et
177  // ----------------------
178  if (jets.size() < 2 || mets.empty()) {
179  preselection = false;
180  }
181 
182  bool write = false;
183 
184  // -----------------------------------------------------
185  //
186  // CompositeCandidates to store the event solution.
187  // This will take one of two forms:
188  // a) lv jj jj Full reconstruction.
189  //
190  // ttbar->
191  // (hadt -> (hadW -> hadp + hadq) + hadb) +
192  // (lept -> (lepW -> lepton + neutrino) + lepb)
193  //
194  // b) lv jj (j) Partial reconstruction, associate
195  // at least 1 jet to the lepton
196  // hemisphere, and at least one jet in
197  // the opposite hemisphere.
198  //
199  // ttbar->
200  // (hadt -> (hadJet1 [+ hadJet2] ) ) +
201  // (lept -> (lepW -> lepton + neutrino) + lepJet1 )
202  //
203  // There will also be two subcategories of (b) that
204  // will correspond to physics cases:
205  //
206  // b1) Lepton is isolated: Moderate ttbar mass.
207  // b2) Lepton is nonisolated: High ttbar mass.
208  //
209  // -----------------------------------------------------
210  reco::CompositeCandidate ttbar("ttbar");
211  AddFourMomenta addFourMomenta;
212 
213  // Main decisions after preselection
214  if (preselection) {
215  if (debug)
216  cout << "Preselection is satisfied" << endl;
217 
218  if (debug)
219  cout << "Jets.size() = " << jets.size() << endl;
220 
221  // This will be modified for the z solution, so make a copy
222  pat::MET neutrino(mets[0]);
223 
224  // 1. First examine the low mass case with 4 jets and widely separated
225  // products. We take out the TtSemiLeptonicEvent from the TQAF and
226  // form the ttbar invariant mass.
227  if (jets.size() >= 4) {
228  if (debug)
229  cout << "Getting ttbar semileptonic solution" << endl;
230 
231  // get the ttbar semileptonic event solution if there are more than 3 jets
233  iEvent.getByToken(solToken_, eSol);
234 
235  // Have solution, continue
236  if (eSol.isValid()) {
237  if (debug)
238  cout << "Got a nonzero size solution vector" << endl;
239  // Just set the ttbar solution to the best ttbar solution from
240  // TtSemiEvtSolutionMaker
242  write = true;
243  }
244  // No ttbar solution with 4 jets, something is weird, print a warning
245  else {
246  edm::LogWarning("DataNotFound") << "BoostedTopProducer: Cannot find TtSemiEvtSolution\n";
247  }
248  }
249  // 2. With 2 or 3 jets, we decide based on the separation between
250  // the lepton and the closest jet in that hemisphere whether to
251  // consider it "moderate" or "high" mass.
252  else if (jets.size() == 2 || jets.size() == 3) {
253  // ------------------------------------------------------------------
254  // First create a leptonic W candidate
255  // ------------------------------------------------------------------
256  reco::CompositeCandidate lepW("lepW");
257 
258  if (isMuon) {
259  if (debug)
260  cout << "Adding muon as daughter" << endl;
261  lepW.addDaughter(*muon, "muon");
262  } else {
263  if (debug)
264  cout << "Adding electron as daughter" << endl;
265  lepW.addDaughter(*electron, "electron");
266  }
267  if (debug)
268  cout << "Adding neutrino as daughter" << endl;
269  lepW.addDaughter(neutrino, "neutrino");
270  addFourMomenta.set(lepW);
271 
272  //bool nuzHasComplex = false;
273  METzCalculator zcalculator;
274 
275  zcalculator.SetMET(neutrino);
276  if (isMuon)
277  zcalculator.SetMuon(*muon);
278  else
279  zcalculator.SetMuon(*electron); // This name is misleading, should be setLepton
280  double neutrinoPz = zcalculator.Calculate(1); // closest to the lepton Pz
281  //if (zcalculator.IsComplex()) nuzHasComplex = true;
282  // Set the neutrino pz
283  neutrino.setPz(neutrinoPz);
284 
285  if (debug)
286  cout << "Set neutrino pz to " << neutrinoPz << endl;
287 
288  // ------------------------------------------------------------------
289  // Next ensure that there is a jet within the hemisphere of the
290  // leptonic W, and one in the opposite hemisphere
291  // ------------------------------------------------------------------
292  reco::CompositeCandidate hadt("hadt");
293  reco::CompositeCandidate lept("lept");
294  if (debug)
295  cout << "Adding lepW as daughter" << endl;
296  lept.addDaughter(lepW, "lepW");
297 
298  std::string hadName("hadJet");
299  std::string lepName("lepJet");
300 
301  // Get the W momentum
302  TLorentzVector p4_W(lepW.px(), lepW.py(), lepW.pz(), lepW.energy());
303 
304  // Loop over the jets
305  std::vector<pat::Jet>::const_iterator jetit = jets.begin(), jetend = jets.end();
306  unsigned long ii = 1; // Count by 1 for naming histograms
307  for (; jetit != jetend; ++jetit, ++ii) {
308  // Get this jet's momentum
309  TLorentzVector p4_jet(jetit->px(), jetit->py(), jetit->pz(), jetit->energy());
310 
311  // Calculate psi (like DeltaR, only more invariant under Rapidity)
312  double psi = Psi(p4_W, p4_jet, mTop_);
313 
314  // Get jets that are in the leptonic hemisphere
315  if (psi < TMath::Pi()) {
316  // Add this jet to the leptonic top
317  std::stringstream s;
318  s << lepName << ii;
319  if (debug)
320  cout << "Adding daughter " << s.str() << endl;
321  lept.addDaughter(*jetit, s.str());
322  }
323  // Get jets that are in the hadronic hemisphere
324  if (psi > TMath::Pi()) {
325  // Add this jet to the hadronic top. We don't
326  // make any W hypotheses in this case, since
327  // we cannot determine which of the three
328  // jets are merged.
329  std::stringstream s;
330  s << hadName << ii;
331  if (debug)
332  cout << "Adding daughter " << s.str() << endl;
333  hadt.addDaughter(*jetit, s.str());
334  }
335  } // end loop over jets
336 
337  addFourMomenta.set(lept);
338  addFourMomenta.set(hadt);
339 
340  bool lepWHasJet = lept.numberOfDaughters() >= 2; // W and >= 1 jet
341  bool hadWHasJet = hadt.numberOfDaughters() >= 1; // >= 1 jet
342  if (lepWHasJet && hadWHasJet) {
343  if (debug)
344  cout << "Adding daughters lept and hadt" << endl;
345  ttbar.addDaughter(lept, "lept");
346  ttbar.addDaughter(hadt, "hadt");
347  addFourMomenta.set(ttbar);
348  write = true;
349  } // end of hadronic jet and leptonic jet
350 
351  } // end if there are 2 or 3 jets
352 
353  } // end if preselection is satisfied
354 
355  // Write the solution to the event record
356  std::vector<reco::CompositeCandidate> ttbarList;
357  if (write) {
358  if (debug)
359  cout << "Writing out" << endl;
360  ttbarList.push_back(ttbar);
361  }
362  std::unique_ptr<std::vector<reco::CompositeCandidate> > pTtbar(new std::vector<reco::CompositeCandidate>(ttbarList));
363  iEvent.put(std::move(pTtbar));
364 }
365 
366 double BoostedTopProducer::Psi(const TLorentzVector& p1, const TLorentzVector& p2, double mass) {
367  TLorentzVector ptot = p1 + p2;
368  Double_t theta1 = TMath::ACos((p1.Vect().Dot(ptot.Vect())) / (p1.P() * ptot.P()));
369  Double_t theta2 = TMath::ACos((p2.Vect().Dot(ptot.Vect())) / (p2.P() * ptot.P()));
370  //Double_t sign = 1.;
371  //if ( (theta1+theta2) > (TMath::Pi()/2) ) sign = -1.;
372  double th1th2 = theta1 + theta2;
373  double psi = (p1.P() + p2.P()) * TMath::Abs(TMath::Sin(th1th2)) / (2. * mass);
374  if (th1th2 > (TMath::Pi() / 2))
375  psi = (p1.P() + p2.P()) * (1. + TMath::Abs(TMath::Cos(th1th2))) / (2. * mass);
376 
377  return psi;
378 }
379 
380 //define this as a plug-in
PDWG_BPHSkim_cff.muons
muons
Definition: PDWG_BPHSkim_cff.py:47
BoostedTopProducer::solToken_
edm::EDGetTokenT< TtSemiLeptonicEvent > solToken_
Definition: BoostedTopProducer.h:95
muon
Definition: MuonCocktails.h:17
BoostedTopProducer.h
BoostedTopProducer::BoostedTopProducer
BoostedTopProducer(const edm::ParameterSet &)
Definition: BoostedTopProducer.cc:23
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
Electron
Definition: Electron.py:1
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
BoostedTopProducer::~BoostedTopProducer
~BoostedTopProducer() override
Definition: BoostedTopProducer.cc:35
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
BoostedTopProducer::eleToken_
edm::EDGetTokenT< std::vector< pat::Electron > > eleToken_
Definition: BoostedTopProducer.h:91
singleTopDQM_cfi.mets
mets
Definition: singleTopDQM_cfi.py:43
BoostedTopProducer::caloIsoCut_
double caloIsoCut_
Definition: BoostedTopProducer.h:98
BoostedTopProducer::muoToken_
edm::EDGetTokenT< std::vector< pat::Muon > > muoToken_
Definition: BoostedTopProducer.h:92
edm::Handle
Definition: AssociativeIterator.h:50
Muon
Definition: Muon.py:1
reco::LeafCandidate::setPz
void setPz(double pz) final
Definition: LeafCandidate.h:163
TtEvent::eventHypo
const reco::CompositeCandidate & eventHypo(const HypoClassKey &key, const unsigned &cmb=0) const
Definition: TtEvent.h:63
TtSemiLeptonicEvent
Class derived from the TtEvent for the semileptonic decay channel.
Definition: TtSemiLeptonicEvent.h:24
alignCSCRings.s
s
Definition: alignCSCRings.py:92
debug
#define debug
Definition: HDRShower.cc:19
BoostedTopProducer::Psi
double Psi(const TLorentzVector &p1, const TLorentzVector &p2, double mass)
Definition: BoostedTopProducer.cc:366
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Jet
Definition: Jet.py:1
Abs
T Abs(T a)
Definition: MathUtil.h:49
BoostedTopProducer
Definition: BoostedTopProducer.h:80
reco::LeafCandidate::py
double py() const final
y coordinate of momentum vector
Definition: LeafCandidate.h:142
BoostedTopProducer::jetToken_
edm::EDGetTokenT< std::vector< pat::Jet > > jetToken_
Definition: BoostedTopProducer.h:93
singleTopDQM_cfi.preselection
preselection
Definition: singleTopDQM_cfi.py:78
metsig::electron
Definition: SignAlgoResolutions.h:48
p2
double p2[4]
Definition: TauolaWrapper.h:90
TtEvent::kMVADisc
Definition: TtEvent.h:35
reco::CompositeCandidate::addDaughter
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
Definition: CompositeCandidate.cc:108
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::LogWarning
Definition: MessageLogger.h:141
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
METzCalculator::Calculate
double Calculate(int type=0)
member functions
Definition: METzCalculator.cc:12
AddFourMomenta.h
iEvent
int iEvent
Definition: GenABIO.cc:224
pat::MET
Analysis-level MET class.
Definition: MET.h:40
p1
double p1[4]
Definition: TauolaWrapper.h:89
edm::EventSetup
Definition: EventSetup.h:57
pat
Definition: HeavyIon.h:7
writeEcalDQMStatus.write
write
Definition: writeEcalDQMStatus.py:48
reco::CompositeCandidate::numberOfDaughters
size_type numberOfDaughters() const override
number of daughters
Definition: CompositeCandidate.cc:41
psi
std::map< std::string, int, std::less< std::string > > psi
Definition: CountProcessesAction.h:15
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
BoostedTopProducer::mTop_
double mTop_
Definition: BoostedTopProducer.h:99
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
BoostedTopProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: BoostedTopProducer.cc:42
BoostedTopProducer::metToken_
edm::EDGetTokenT< std::vector< pat::MET > > metToken_
Definition: BoostedTopProducer.h:94
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
AddFourMomenta::set
void set(reco::Candidate &c) const
set up a candidate
Definition: AddFourMomenta.cc:6
MET
reco::LeafCandidate::energy
double energy() const final
energy
Definition: LeafCandidate.h:125
Pi
const double Pi
Definition: CosmicMuonParameters.h:18
AddFourMomenta
Definition: AddFourMomenta.h:18
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
edm::Event
Definition: Event.h:73
reco::LeafCandidate::px
double px() const final
x coordinate of momentum vector
Definition: LeafCandidate.h:140
METzCalculator::SetMET
void SetMET(const pat::MET &MET)
Set MET.
Definition: METzCalculator.h:29
reco::LeafCandidate::pz
double pz() const final
z coordinate of momentum vector
Definition: LeafCandidate.h:144
METzCalculator
Definition: METzCalculator.h:21
cuy.ii
ii
Definition: cuy.py:590
reco::CompositeCandidate
Definition: CompositeCandidate.h:21
METzCalculator::SetMuon
void SetMuon(const pat::Particle &lepton)
Set Muon.
Definition: METzCalculator.h:35
reco::isMuon
bool isMuon(const Candidate &part)
Definition: pdgIdUtils.h:9