CMS 3D CMS Logo

TTbar_Kinematics.cc
Go to the documentation of this file.
6 
7 using namespace edm;
9  : hepmcCollection_(iConfig.getParameter<edm::InputTag>("hepmcCollection")),
10  genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag")) {
11  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
12  genEventInfoProductTagToken_ = consumes<GenEventInfoProduct>(genEventInfoProductTag_);
13 }
14 
16 
17 //
18 // member functions
19 //
20 
21 // ------------ method called for each event ------------
23  // --- the MC weights ---
25  iEvent.getByToken(genEventInfoProductTagToken_, evt_info);
26  if (!evt_info.isValid())
27  return;
28  weight = evt_info->weight();
29 
32  iEvent.getByToken(hepmcCollectionToken_, evt);
33 
34  //Get EVENT
35  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
36 
38  bool top(false), antitop(false), antibottom(false), bottom(false), Wplus(false), Wmin(false);
39  for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
40  iter != myGenEvent->particles_end();
41  iter++) {
42  if ((*iter)->pdg_id() == PdtPdgMini::t || (*iter)->pdg_id() == PdtPdgMini::anti_t) {
43  if ((*iter)->end_vertex()) {
44  HepMC::GenVertex::particle_iterator des;
45  for (des = (*iter)->end_vertex()->particles_begin(HepMC::children);
46  des != (*iter)->end_vertex()->particles_end(HepMC::children);
47  ++des) {
48  if ((*des)->pdg_id() == PdtPdgMini::b) {
49  tlv_Bottom.SetPxPyPzE(
50  (*des)->momentum().px(), (*des)->momentum().py(), (*des)->momentum().pz(), (*des)->momentum().e());
51  bottom = true;
52  }
53  if ((*des)->pdg_id() == PdtPdgMini::anti_b) {
54  antibottom = true;
55  tlv_BottomBar.SetPxPyPzE(
56  (*des)->momentum().px(), (*des)->momentum().py(), (*des)->momentum().pz(), (*des)->momentum().e());
57  }
58  if ((*des)->pdg_id() == PdtPdgMini::W_plus) {
59  tlv_TopBar.SetPxPyPzE(
60  (*iter)->momentum().px(), (*iter)->momentum().py(), (*iter)->momentum().pz(), (*iter)->momentum().e());
61  antitop = true;
62  tlv_Wplus.SetPxPyPzE(
63  (*des)->momentum().px(), (*des)->momentum().py(), (*des)->momentum().pz(), (*des)->momentum().e());
64  Wplus = true;
65  }
66  if ((*des)->pdg_id() == PdtPdgMini::W_minus) {
67  tlv_Top.SetPxPyPzE(
68  (*iter)->momentum().px(), (*iter)->momentum().py(), (*iter)->momentum().pz(), (*iter)->momentum().e());
69  top = true;
70  tlv_Wmin.SetPxPyPzE(
71  (*des)->momentum().px(), (*des)->momentum().py(), (*des)->momentum().pz(), (*des)->momentum().e());
72  Wmin = true;
73  }
74  }
75  }
76  }
77  }
78 
80 
81  //---topquarkquantities---
82  nEvt->Fill(0.5, weight);
83  if (top && antitop) {
84  hTopPt->Fill(tlv_Top.Pt(), weight);
85  hTopPt->Fill(tlv_TopBar.Pt(), weight);
86 
87  hTopY->Fill(tlv_Top.Rapidity(), weight);
88  hTopY->Fill(tlv_TopBar.Rapidity(), weight);
89 
90  hTopMass->Fill(tlv_Top.M(), weight);
92 
93  //---ttbarpairquantities---
94  hTTbarPt->Fill(tlv_TTbar.Pt(), weight);
95  hTTbarY->Fill(tlv_TTbar.Rapidity(), weight);
97  }
98  if (bottom && antibottom) {
101 
102  hBottomEta->Fill(tlv_Bottom.Eta(), weight);
104 
105  //hBottomY->Fill(math::XYZTLorentzVector(bottom->momentum()).Rapidity(),weight);
106  //hBottomY->Fill(math::XYZTLorentzVector(antibottom->momentum()).Rapidity(),weight);
107 
108  hBottomY->Fill(tlv_Bottom.Rapidity(), weight);
109  hBottomY->Fill(tlv_BottomBar.Rapidity(), weight);
110 
113 
116 
119 
122 
125 
128 
129  hBottomYPz->Fill(tlv_Bottom.Rapidity(), tlv_Bottom.Pz(), weight);
130  hBottomYPz->Fill(tlv_BottomBar.Rapidity(), tlv_BottomBar.Pz(), weight);
131 
134 
137 
138  hBottomMassY->Fill(tlv_Bottom.M(), tlv_Bottom.Rapidity(), weight);
140 
141  hBottomMassDeltaY->Fill(tlv_Bottom.M(), tlv_Bottom.Eta() - tlv_Bottom.Rapidity(), weight);
143  }
144  if (Wplus && Wmin) {
145  hWplusPz->Fill(tlv_Wplus.Pz(), weight);
146  hWminPz->Fill(tlv_Wmin.Pz(), weight);
147  }
148 }
149 
150 // ------------ method called once each job just before starting event loop ------------
152  DQMHelper dqm(&i);
153  i.setCurrentFolder("Generator/TTbar");
154 
155  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "bins", "Number of Events");
156 
157  hTopPt = dqm.book1dHisto(
158  "TTbar_TopPt", "t quark transverse momentum", 1000, 0., 1000., "P_{t}^{t quark} (GeV)", "Number of Events");
159  hTopY = dqm.book1dHisto("TTbar_TopY", "t quark rapidity", 200, -5., 5., "Y_{t quark}", "Number of Events");
160  hTopMass = dqm.book1dHisto("TTbar_TopMass", "t quark mass", 500, 0., 500., "M_{t quark} (GeV)", "Number of Events");
161 
162  hTTbarPt = dqm.book1dHisto(
163  "TTbar_TTbarPt", "tt pair transverse momentum", 1000, 0., 1000., "P_{t}^{tt pair} (GeV)", "Number of Events");
164  hTTbarY = dqm.book1dHisto("TTbar_TTbarY", "tt pair rapidity", 200, -5., 5., "Y_{tt pair}", "Number of Events");
165  hTTbarMass =
166  dqm.book1dHisto("TTbar_TTbarMass", "tt pair mass", 1000, 0., 1000., "M_{tt pair} (GeV)", "Number of Events");
167 
168  hBottomPt = dqm.book1dHisto(
169  "TTbar_BottomPt", "b quark transverse momentum", 1000, 0., 1000., "P_{t}^{b quark} (GeV)", "Number of Events");
170  hBottomEta = dqm.book1dHisto(
171  "TTbar_BottomEta", "b quark pseudo-rapidity", 200, -5., 5., "#eta_{b quark} (GeV)", "Number of Events");
172  hBottomY =
173  dqm.book1dHisto("TTbar_BottomY", "b quark rapidity", 200, -5., 5., "M_{b quark} (GeV)", "Number of Events");
174  hBottomPz = dqm.book1dHisto(
175  "TTbar_BottomPz", "b quark longitudinal momentum", 200, -100., 100., "P_{z}^{b quark} (GeV)", "Number of Events");
176  hBottomE =
177  dqm.book1dHisto("TTbar_BottomE", "b quark energy", 1000, 0., 1000., "E_{b quark} (GeV)", "Number of Events");
178  hBottomMass =
179  dqm.book1dHisto("TTbar_BottomMass", "b quark mass", 50, 0., 5., "M_{b quark} (GeV)", "Number of Events");
180 
181  hBottomPtPz = dqm.book2dHisto("TTbar_BottomPtPz",
182  "b quark longitudinal vs transverse momentum",
183  1000,
184  0.,
185  1000.,
186  200,
187  -100.,
188  100.,
189  "P_{z}^{b quark} (GeV)",
190  "P_{t}^{b quark} (GeV)");
191  hBottomEtaPz = dqm.book2dHisto("TTbar_BottomEtaPz",
192  "b quark longitudinal momentum vs pseudorapidity",
193  200,
194  -5.,
195  5.,
196  200,
197  -100.,
198  100.,
199  "#eta_{b quark}",
200  "P_{z}^{b quark} (GeV)");
201  hBottomEtaPt = dqm.book2dHisto("TTbar_BottomEtaPt",
202  " quark transveral momentum vs pseudorapidity",
203  200,
204  -5.,
205  5.,
206  1000,
207  0.,
208  1000.,
209  "#eta_{b quark}",
210  "P_{t}^{b quark} (GeV)");
211  hBottomYPz = dqm.book2dHisto("TTbar_BottomYPz",
212  "b quark longitudinal momentum vs rapidity",
213  200,
214  -5.,
215  5.,
216  200,
217  -100.,
218  100.,
219  "Y_{b quark}",
220  "P_{z}^{b quark} (GeV)");
221  hBottomMassPz = dqm.book2dHisto("TTbar_BottomMassPz",
222  "b quark longitudinal momentum vs mass",
223  50,
224  0.,
225  5.,
226  200,
227  -100.,
228  100.,
229  "M_{b quark} (GeV)",
230  "P_{z}^{b quark} (GeV)");
231  hBottomMassEta = dqm.book2dHisto("TTbar_BottomMassEta",
232  "b quark pseudorapidity vs mass",
233  50,
234  0.,
235  5.,
236  200,
237  -5.,
238  5.,
239  "M_{b quark} (GeV)",
240  "#eta_{b quark}");
241  hBottomMassY = dqm.book2dHisto(
242  "TTbar_BottomMassY", "b quark rapidity vs mass", 50, 0., 5., 200, -5., 5., "M_{b quark} (GeV)", "Y_{b quark}");
243  hBottomMassDeltaY = dqm.book2dHisto("TTbar_BottomMassDeltaY",
244  "b quark pseudorapidity - rapidity vs mass",
245  50,
246  0.,
247  50.,
248  2000,
249  -5.,
250  5.,
251  "M_{b quark} (GeV)",
252  "Y_{b quark}");
253 
254  hWplusPz = dqm.book1dHisto(
255  "TTbar_WplusPz", "W+ boson longitudinal momentum", 200, -100., 100., "P_{z}^{W+} (GeV)", "Number of Events");
256  hWminPz = dqm.book1dHisto(
257  "TTbar_WminPz", "W- boson longitudinal momentum", 200, -100., 100., "P_{z}^{W-} (GeV)", "Number of Events");
258 }
PdtPdgMini.h
TTbar_Kinematics::TTbar_Kinematics
TTbar_Kinematics(const edm::ParameterSet &)
Definition: TTbar_Kinematics.cc:8
mps_fire.i
i
Definition: mps_fire.py:428
TTbar_Kinematics::hBottomMass
MonitorElement * hBottomMass
Definition: TTbar_Kinematics.h:92
TTbar_Kinematics::hBottomYPz
MonitorElement * hBottomYPz
Definition: TTbar_Kinematics.h:100
PdtPdgMini::W_plus
Definition: PdtPdgMini.h:42
edm::Run
Definition: Run.h:45
TTbar_Kinematics::tlv_Bottom
TLorentzVector tlv_Bottom
Definition: TTbar_Kinematics.h:71
TTbar_Kinematics::tlv_Wplus
TLorentzVector tlv_Wplus
Definition: TTbar_Kinematics.h:73
edm
HLT enums.
Definition: AlignableModifier.h:19
class-composition.children
children
Definition: class-composition.py:98
TTbar_Kinematics::hBottomY
MonitorElement * hBottomY
Definition: TTbar_Kinematics.h:89
TTbar_Kinematics::bookHistograms
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Definition: TTbar_Kinematics.cc:151
PdtPdgMini::W_minus
Definition: PdtPdgMini.h:43
TTbar_Kinematics::hBottomEtaPt
MonitorElement * hBottomEtaPt
Definition: TTbar_Kinematics.h:99
DQMStore.h
TTbar_Kinematics::hBottomEta
MonitorElement * hBottomEta
Definition: TTbar_Kinematics.h:88
TTbar_Kinematics::~TTbar_Kinematics
~TTbar_Kinematics() override
Definition: TTbar_Kinematics.cc:15
TTbar_Kinematics::genEventInfoProductTagToken_
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoProductTagToken_
Definition: TTbar_Kinematics.h:106
edm::Handle
Definition: AssociativeIterator.h:50
TTbar_Kinematics::hBottomPtPz
MonitorElement * hBottomPtPz
Definition: TTbar_Kinematics.h:97
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
TTbar_Kinematics::hTopMass
MonitorElement * hTopMass
Definition: TTbar_Kinematics.h:81
GenParticle.h
DQMHelper.h
TTbar_Kinematics::hBottomMassY
MonitorElement * hBottomMassY
Definition: TTbar_Kinematics.h:103
PdtPdgMini::anti_b
Definition: PdtPdgMini.h:16
TTbar_Kinematics::hTTbarMass
MonitorElement * hTTbarMass
Definition: TTbar_Kinematics.h:85
TTbar_Kinematics.h
PdtPdgMini::anti_t
Definition: PdtPdgMini.h:18
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition: MonitorElement.h:290
TTbar_Kinematics::tlv_TopBar
TLorentzVector tlv_TopBar
Definition: TTbar_Kinematics.h:70
TTbar_Kinematics::weight
double weight
Definition: TTbar_Kinematics.h:67
TTbar_Kinematics::hBottomPz
MonitorElement * hBottomPz
Definition: TTbar_Kinematics.h:90
TTbar_Kinematics::hTopY
MonitorElement * hTopY
Definition: TTbar_Kinematics.h:80
PdtPdgMini::b
Definition: PdtPdgMini.h:15
edm::ParameterSet
Definition: ParameterSet.h:47
TTbar_Kinematics::tlv_BottomBar
TLorentzVector tlv_BottomBar
Definition: TTbar_Kinematics.h:72
TTbar_Kinematics::hBottomMassPz
MonitorElement * hBottomMassPz
Definition: TTbar_Kinematics.h:101
TTbar_Kinematics::hBottomPt
MonitorElement * hBottomPt
Definition: TTbar_Kinematics.h:87
TTbar_Kinematics::genEventInfoProductTag_
edm::InputTag genEventInfoProductTag_
Definition: TTbar_Kinematics.h:65
PdtPdgMini::t
Definition: PdtPdgMini.h:17
iEvent
int iEvent
Definition: GenABIO.cc:224
TTbar_Kinematics::hTTbarY
MonitorElement * hTTbarY
Definition: TTbar_Kinematics.h:84
edm::EventSetup
Definition: EventSetup.h:58
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
TTbar_Kinematics::tlv_TTbar
TLorentzVector tlv_TTbar
Definition: TTbar_Kinematics.h:76
TTbar_Kinematics::hBottomMassDeltaY
MonitorElement * hBottomMassDeltaY
Definition: TTbar_Kinematics.h:104
TTbar_Kinematics::hWminPz
MonitorElement * hWminPz
Definition: TTbar_Kinematics.h:95
TTbar_Kinematics::hBottomE
MonitorElement * hBottomE
Definition: TTbar_Kinematics.h:91
alignCSCRings.r
r
Definition: alignCSCRings.py:93
DQMHelper
Definition: DQMHelper.h:15
TTbar_Kinematics::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: TTbar_Kinematics.cc:22
TTbar_Kinematics::tlv_Wmin
TLorentzVector tlv_Wmin
Definition: TTbar_Kinematics.h:74
TTbar_Kinematics::hTTbarPt
MonitorElement * hTTbarPt
Definition: TTbar_Kinematics.h:83
TTbar_Kinematics::tlv_Top
TLorentzVector tlv_Top
Definition: TTbar_Kinematics.h:69
TTbar_Kinematics::hBottomMassEta
MonitorElement * hBottomMassEta
Definition: TTbar_Kinematics.h:102
TTbar_Kinematics::hTopPt
MonitorElement * hTopPt
Definition: TTbar_Kinematics.h:79
dqm::implementation::IBooker
Definition: DQMStore.h:43
dqm
Definition: DQMStore.h:18
TTbar_Kinematics::hBottomEtaPz
MonitorElement * hBottomEtaPz
Definition: TTbar_Kinematics.h:98
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
edm::Event
Definition: Event.h:73
TTbar_Kinematics::nEvt
MonitorElement * nEvt
Definition: TTbar_Kinematics.h:78
TTbar_Kinematics::hWplusPz
MonitorElement * hWplusPz
Definition: TTbar_Kinematics.h:94
TTbar_Kinematics::hepmcCollectionToken_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
Definition: TTbar_Kinematics.h:107
edm::InputTag
Definition: InputTag.h:15
TTbar_Kinematics::hepmcCollection_
edm::InputTag hepmcCollection_
Definition: TTbar_Kinematics.h:64
weight
Definition: weight.py:1
GenEventInfoProduct::weight
double weight() const
Definition: GenEventInfoProduct.h:35
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37