CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TTbar_Kinematics.cc
Go to the documentation of this file.
5 
6 
7 using namespace edm;
9  hepmcCollection_(iConfig.getParameter<edm::InputTag>("hepmcCollection"))
10  ,genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag"))
11 {
12  dbe = 0;
14 
15  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
16  genEventInfoProductTagToken_=consumes<GenEventInfoProduct>(genEventInfoProductTag_);
17 
18 }
19 
20 
22 {
23 
24  // do anything here that needs to be done at desctruction time
25  // (e.g. close files, deallocate resources etc.)
26 
27 }
28 
29 
30 //
31 // member functions
32 //
33 
34 // ------------ method called for each event ------------
35 void
37 {
38 
39  // --- the MC weights ---
41  iEvent.getByToken(genEventInfoProductTagToken_, evt_info);
42  if(!evt_info.isValid()) return;
43  weight = evt_info->weight() ;
44 
45 
46 
49  iEvent.getByToken(hepmcCollectionToken_, evt);
50 
51  //Get EVENT
52  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
53 
55  bool top(false), antitop(false), antibottom(false), bottom(false), Wplus(false), Wmin(false);
56  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); iter++) {
57  if((*iter)->pdg_id()==PdtPdgMini::t || (*iter)->pdg_id()==PdtPdgMini::anti_t){
58  if( (*iter)->end_vertex()){
59  HepMC::GenVertex::particle_iterator des;
60  for(des = (*iter)->end_vertex()->particles_begin(HepMC::children);des!= (*iter)->end_vertex()->particles_end(HepMC::children);++des ){
61  if((*des)->pdg_id()==PdtPdgMini::b){
62  tlv_Bottom.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
63  bottom=true;
64  }
65  if((*des)->pdg_id()==PdtPdgMini::anti_b){
66  antibottom=true;
67  tlv_BottomBar.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
68  }
69  if((*des)->pdg_id()==PdtPdgMini::W_plus){
70  tlv_TopBar.SetPxPyPzE((*iter)->momentum().px(),(*iter)->momentum().py(),(*iter)->momentum().pz(),(*iter)->momentum().e()); antitop=true;
71  tlv_Wplus.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e()); Wplus=true;
72  }
73  if((*des)->pdg_id()==PdtPdgMini::W_minus){
74  tlv_Top.SetPxPyPzE((*iter)->momentum().px(),(*iter)->momentum().py(),(*iter)->momentum().pz(),(*iter)->momentum().e()); top=true;
75  tlv_Wmin.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e()); Wmin=true;
76  }
77  }
78  }
79  }
80  }
81 
82  tlv_TTbar = tlv_Top + tlv_TopBar ;
83 
84  //---topquarkquantities---
85  nEvt->Fill(0.5,weight);
86  if(top && antitop){
87  hTopPt->Fill(tlv_Top.Pt(),weight);
88  hTopPt->Fill(tlv_TopBar.Pt(),weight);
89 
90  hTopY->Fill(tlv_Top.Rapidity(),weight);
91  hTopY->Fill(tlv_TopBar.Rapidity(),weight);
92 
93  hTopMass->Fill(tlv_Top.M(),weight);
94  hTopMass->Fill(tlv_TopBar.M(),weight);
95 
96  //---ttbarpairquantities---
97  hTTbarPt->Fill(tlv_TTbar.Pt(),weight);
98  hTTbarY->Fill(tlv_TTbar.Rapidity(),weight);
99  hTTbarMass->Fill(tlv_TTbar.M(),weight);
100  }
101  if(bottom && antibottom){
102  hBottomPt->Fill(tlv_Bottom.Pt(),weight);
103  hBottomPt->Fill(tlv_BottomBar.Pt(),weight);
104 
105  hBottomEta->Fill(tlv_Bottom.Eta(),weight);
106  hBottomEta->Fill(tlv_BottomBar.Eta(),weight);
107 
108  //hBottomY->Fill(math::XYZTLorentzVector(bottom->momentum()).Rapidity(),weight);
109  //hBottomY->Fill(math::XYZTLorentzVector(antibottom->momentum()).Rapidity(),weight);
110 
111  hBottomY->Fill(tlv_Bottom.Rapidity(),weight);
112  hBottomY->Fill(tlv_BottomBar.Rapidity(),weight);
113 
114  hBottomPz->Fill(tlv_Bottom.Pz(),weight);
115  hBottomPz->Fill(tlv_BottomBar.Pz(),weight);
116 
117  hBottomE->Fill(tlv_Bottom.E(),weight);
118  hBottomE->Fill(tlv_BottomBar.E(),weight);
119 
120  hBottomMass->Fill(tlv_Bottom.M(),weight);
121  hBottomMass->Fill(tlv_BottomBar.M(),weight);
122 
123  hBottomPtPz->Fill(tlv_Bottom.Pt(),tlv_Bottom.Pz(),weight);
124  hBottomPtPz->Fill(tlv_BottomBar.Pt(),tlv_BottomBar.Pz(),weight);
125 
126  hBottomEtaPz->Fill(tlv_Bottom.Eta(),tlv_Bottom.Pz(),weight);
127  hBottomEtaPz->Fill(tlv_BottomBar.Eta(),tlv_BottomBar.Pz(),weight);
128 
129  hBottomEtaPt->Fill(tlv_Bottom.Eta(),tlv_Bottom.Pt(),weight);
130  hBottomEtaPt->Fill(tlv_BottomBar.Eta(),tlv_BottomBar.Pt(),weight);
131 
132  hBottomYPz->Fill(tlv_Bottom.Rapidity(),tlv_Bottom.Pz(),weight);
133  hBottomYPz->Fill(tlv_BottomBar.Rapidity(),tlv_BottomBar.Pz(),weight);
134 
135  hBottomMassPz->Fill(tlv_Bottom.M(),tlv_Bottom.Pz(),weight);
136  hBottomMassPz->Fill(tlv_BottomBar.M(),tlv_BottomBar.Pz(),weight);
137 
138  hBottomMassEta->Fill(tlv_Bottom.M(),tlv_Bottom.Eta(),weight);
139  hBottomMassEta->Fill(tlv_BottomBar.M(),tlv_BottomBar.Eta(),weight);
140 
141  hBottomMassY->Fill(tlv_Bottom.M(),tlv_Bottom.Rapidity(),weight);
142  hBottomMassY->Fill(tlv_BottomBar.M(),tlv_BottomBar.Rapidity(),weight);
143 
144  hBottomMassDeltaY->Fill(tlv_Bottom.M(),tlv_Bottom.Eta()-tlv_Bottom.Rapidity(),weight);
145  hBottomMassDeltaY->Fill(tlv_BottomBar.M(),tlv_BottomBar.Eta()-tlv_BottomBar.Rapidity(),weight);
146  }
147  if(Wplus && Wmin){
148  hWplusPz->Fill(tlv_Wplus.Pz(),weight);
149  hWminPz->Fill(tlv_Wmin.Pz(),weight);
150  }
151 }
152 
153 
154 // ------------ method called once each job just before starting event loop ------------
155 void
157 {
158  if(!dbe) return;
159  dbe->setCurrentFolder("Generator/TTbar");
160 
161  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
162 
163  hTopPt = dbe->book1D("TTbar_TopPt","t quark transverse momentum",1000,0.,1000.); hTopPt->setAxisTitle("t quark transverse momentum");
164  hTopY = dbe->book1D("TTbar_TopY","t quark rapidity",200,-5.,5.); hTopY->setAxisTitle("t quark rapidity");
165  hTopMass = dbe->book1D("TTbar_TopMass","t quark mass",500,0.,500.); hTopMass->setAxisTitle("t quark mass");
166 
167  hTTbarPt = dbe->book1D("TTbar_TTbarPt","tt pair transverse momentum",1000,0.,1000.); hTTbarPt->setAxisTitle("tt pair transverse momentum");
168  hTTbarY = dbe->book1D("TTbar_TTbarY","tt pair rapidity",200,-5.,5.); hTTbarY->setAxisTitle("tt pair rapidity");
169  hTTbarMass = dbe->book1D("TTbar_TTbarMass","tt pair mass",1000,0.,1000.); hTTbarMass->setAxisTitle("tt pair mass");
170 
171  hBottomPt = dbe->book1D("TTbar_BottomPt","b quark transverse momentum",1000,0.,1000.); hBottomPt->setAxisTitle("b quark transverse momentum");
172  hBottomEta = dbe->book1D("TTbar_BottomEta","b quark pseudo-rapidity",200,-5.,5.); hBottomEta->setAxisTitle("b quark pseudo-rapidity");
173  hBottomY = dbe->book1D("TTbar_BottomY","b quark rapidity",200,-5.,5.); hBottomY->setAxisTitle("b quark rapidity");
174  hBottomPz = dbe->book1D("TTbar_BottomPz","b quark longitudinal momentum",200,-100.,100.); hBottomPz->setAxisTitle("b quark longitudinal momentum");
175  hBottomE = dbe->book1D("TTbar_BottomE","b quark energy",1000,0.,1000.); hBottomE->setAxisTitle("b quark energy");
176  hBottomMass = dbe->book1D("TTbar_BottomMass","b quark mass",50,0.,5.); hBottomMass->setAxisTitle("b quark mass");
177 
178  hBottomPtPz = dbe->book2D("TTbar_BottomPtPz","b quark longitudinal vs transverse momentum",1000,0.,1000.,200,-100.,100.); hBottomPtPz->setAxisTitle("P_{z} (GeV)",1);hBottomPtPz->setAxisTitle("P_{t} (GeV)",2);
179  hBottomEtaPz = dbe->book2D("TTbar_BottomEtaPz","b quark longitudinal momentum vs pseudorapidity",200,-5.,5.,200,-100.,100.); hBottomEtaPz->setAxisTitle("#eta",1);hBottomEtaPz->setAxisTitle("P_{z} (GeV)",1);
180  hBottomEtaPt = dbe->book2D("TTbar_BottomEtaPt"," quark transveral momentum vs pseudorapidity",200,-5.,5.,1000,0.,1000.); hBottomEtaPt->setAxisTitle("#eta");hBottomEtaPt->setAxisTitle("P_{t} (GeV)");
181  hBottomYPz = dbe->book2D("TTbar_BottomYPz","b quark longitudinal momentum vs rapidity",200,-5.,5.,200,-100.,100.); hBottomYPz->setAxisTitle("Y");hBottomYPz->setAxisTitle("P_{z} (GeV)");
182  hBottomMassPz = dbe->book2D("TTbar_BottomMassPz","b quark longitudinal momentum vs mass",50,0.,5.,200,-100.,100.); hBottomMassPz->setAxisTitle("M (GeV)");hBottomMassPz->setAxisTitle("P_{z} (GeV)");
183  hBottomMassEta = dbe->book2D("TTbar_BottomMassEta","b quark pseudorapidity vs mass",50,0.,5.,200,-5.,5.); hBottomMassEta->setAxisTitle("M (GeV)");hBottomMassEta->setAxisTitle("#eta");
184  hBottomMassY = dbe->book2D("TTbar_BottomMassY","b quark rapidity vs mass",50,0.,5.,200,-5.,5.); hBottomMassY->setAxisTitle("M (GeV)"); hBottomMassY->setAxisTitle("Y");
185  hBottomMassDeltaY = dbe->book2D("TTbar_BottomMassDeltaY","b quark pseudorapidity - rapidity vs mass",50,0.,50.,2000,-5.,5.); hBottomMassDeltaY->setAxisTitle("M (GeV)");hBottomMassDeltaY->setAxisTitle("Y");
186 
187  hWplusPz = dbe->book1D("TTbar_WplusPz","W+ boson longitudinal momentum",200,-100.,100.); hWplusPz->setAxisTitle("W+ boson longitudinal momentum");
188  hWminPz = dbe->book1D("TTbar_WminPz","W- boson longitudinal momentum",200,-100.,100.); hWminPz->setAxisTitle("W- boson longitudinal momentum");
189 
190 }
191 
192 // ------------ method called once each job just after ending the event loop ------------
193 void
195 {
196 }
197 
198 // ------------ method called when starting to processes a run ------------
199 void
201 {
202 }
203 
204 // ------------ method called when ending the processing of a run ------------
205 void
207 {
208 }
209 
210 // ------------ method called when starting to processes a luminosity block ------------
211 void
213 {
214 }
215 
216 // ------------ method called when ending the processing of a luminosity block ------------
217 void
219 {
220 }
221 
222 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
223 void
225  //The following says we do not know what parameters are allowed so do no validation
226  // Please change this to state exactly what you do use, even if it is no parameters
228  desc.setUnknown();
229  descriptions.addDefault(desc);
230 }
MonitorElement * hBottomYPz
MonitorElement * nEvt
TLorentzVector tlv_Bottom
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoProductTagToken_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:873
MonitorElement * hBottomMass
MonitorElement * hBottomPtPz
virtual void endRun(edm::Run const &, edm::EventSetup const &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
MonitorElement * hTTbarMass
virtual void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * hBottomMassY
TLorentzVector tlv_Wplus
MonitorElement * hBottomEtaPt
virtual void endJob()
void Fill(long long x)
virtual void beginJob()
MonitorElement * hBottomEta
TLorentzVector tlv_TopBar
int iEvent
Definition: GenABIO.cc:243
void addDefault(ParameterSetDescription const &psetDescription)
MonitorElement * hTopMass
MonitorElement * hTTbarY
TLorentzVector tlv_TTbar
MonitorElement * hWminPz
bool isValid() const
Definition: HandleBase.h:76
TLorentzVector tlv_BottomBar
MonitorElement * hBottomMassPz
MonitorElement * hBottomMassEta
edm::InputTag genEventInfoProductTag_
MonitorElement * hBottomPz
MonitorElement * hTopY
TLorentzVector tlv_Wmin
MonitorElement * hTTbarPt
MonitorElement * hBottomMassDeltaY
MonitorElement * hBottomPt
MonitorElement * hTopPt
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
TLorentzVector tlv_Top
MonitorElement * hBottomE
MonitorElement * hWplusPz
DQMStore * dbe
ME&#39;s &quot;container&quot;.
virtual void beginRun(edm::Run const &, edm::EventSetup const &)
edm::InputTag hepmcCollection_
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1001
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * hBottomEtaPz
TTbar_Kinematics(const edm::ParameterSet &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:585
Definition: Run.h:41
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * hBottomY