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 
16 
18 {
19 
20  // do anything here that needs to be done at desctruction time
21  // (e.g. close files, deallocate resources etc.)
22 
23 }
24 
25 
26 //
27 // member functions
28 //
29 
30 // ------------ method called for each event ------------
31 void
33 {
34 
35  // --- the MC weights ---
37  iEvent.getByLabel(genEventInfoProductTag_, evt_info);
38  if(!evt_info.isValid()) return;
39  weight = evt_info->weight() ;
40 
41  /*
42  // --- get TopQuarkAnalysis TtGenEvent
43  Handle<TtGenEvent> genEvt;
44  iEvent.getByLabel("genEvt", genEvt);
45 
46  if(!genEvt.isValid())return;
47  */
48 
51  iEvent.getByLabel(hepmcCollection_, evt);
52 
53  //Get EVENT
54  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
55 
57  bool top(false), antitop(false), antibottom(false), bottom(false), Wplus(false), Wmin(false);
58  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); iter++) {
59  if((*iter)->pdg_id()==PdtPdgMini::t || (*iter)->pdg_id()==PdtPdgMini::anti_t){
60  if( (*iter)->end_vertex()){
61  HepMC::GenVertex::particle_iterator des;
62  for(des = (*iter)->end_vertex()->particles_begin(HepMC::children);des!= (*iter)->end_vertex()->particles_end(HepMC::children);++des ){
63  if((*des)->pdg_id()==PdtPdgMini::b){
64  tlv_Bottom.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
65  bottom=true;
66  }
67  if((*des)->pdg_id()==PdtPdgMini::anti_b){
68  antibottom=true;
69  tlv_BottomBar.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
70  }
71  if((*des)->pdg_id()==PdtPdgMini::W_plus){
72  tlv_TopBar.SetPxPyPzE((*iter)->momentum().px(),(*iter)->momentum().py(),(*iter)->momentum().pz(),(*iter)->momentum().e()); antitop=true;
73  tlv_Wplus.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e()); Wplus=true;
74  }
75  if((*des)->pdg_id()==PdtPdgMini::W_minus){
76  tlv_Top.SetPxPyPzE((*iter)->momentum().px(),(*iter)->momentum().py(),(*iter)->momentum().pz(),(*iter)->momentum().e()); top=true;
77  tlv_Wmin.SetPxPyPzE((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e()); Wmin=true;
78  }
79  }
80  }
81  }
82  }
83 
84  tlv_TTbar = tlv_Top + tlv_TopBar ;
85 
86  //---topquarkquantities---
87  nEvt->Fill(0.5,weight);
88  if(top && antitop){
89  hTopPt->Fill(tlv_Top.Pt(),weight);
90  hTopPt->Fill(tlv_TopBar.Pt(),weight);
91 
92  hTopY->Fill(tlv_Top.Rapidity(),weight);
93  hTopY->Fill(tlv_TopBar.Rapidity(),weight);
94 
95  hTopMass->Fill(tlv_Top.M(),weight);
96  hTopMass->Fill(tlv_TopBar.M(),weight);
97 
98  //---ttbarpairquantities---
99  hTTbarPt->Fill(tlv_TTbar.Pt(),weight);
100  hTTbarY->Fill(tlv_TTbar.Rapidity(),weight);
101  hTTbarMass->Fill(tlv_TTbar.M(),weight);
102  }
103  if(bottom && antibottom){
104  hBottomPt->Fill(tlv_Bottom.Pt(),weight);
105  hBottomPt->Fill(tlv_BottomBar.Pt(),weight);
106 
107  hBottomEta->Fill(tlv_Bottom.Eta(),weight);
108  hBottomEta->Fill(tlv_BottomBar.Eta(),weight);
109 
110  //hBottomY->Fill(math::XYZTLorentzVector(bottom->momentum()).Rapidity(),weight);
111  //hBottomY->Fill(math::XYZTLorentzVector(antibottom->momentum()).Rapidity(),weight);
112 
113  hBottomY->Fill(tlv_Bottom.Rapidity(),weight);
114  hBottomY->Fill(tlv_BottomBar.Rapidity(),weight);
115 
116  hBottomPz->Fill(tlv_Bottom.Pz(),weight);
117  hBottomPz->Fill(tlv_BottomBar.Pz(),weight);
118 
119  hBottomE->Fill(tlv_Bottom.E(),weight);
120  hBottomE->Fill(tlv_BottomBar.E(),weight);
121 
122  hBottomMass->Fill(tlv_Bottom.M(),weight);
123  hBottomMass->Fill(tlv_BottomBar.M(),weight);
124 
125  hBottomPtPz->Fill(tlv_Bottom.Pt(),tlv_Bottom.Pz(),weight);
126  hBottomPtPz->Fill(tlv_BottomBar.Pt(),tlv_BottomBar.Pz(),weight);
127 
128  hBottomEtaPz->Fill(tlv_Bottom.Eta(),tlv_Bottom.Pz(),weight);
129  hBottomEtaPz->Fill(tlv_BottomBar.Eta(),tlv_BottomBar.Pz(),weight);
130 
131  hBottomEtaPt->Fill(tlv_Bottom.Eta(),tlv_Bottom.Pt(),weight);
132  hBottomEtaPt->Fill(tlv_BottomBar.Eta(),tlv_BottomBar.Pt(),weight);
133 
134  hBottomYPz->Fill(tlv_Bottom.Rapidity(),tlv_Bottom.Pz(),weight);
135  hBottomYPz->Fill(tlv_BottomBar.Rapidity(),tlv_BottomBar.Pz(),weight);
136 
137  hBottomMassPz->Fill(tlv_Bottom.M(),tlv_Bottom.Pz(),weight);
138  hBottomMassPz->Fill(tlv_BottomBar.M(),tlv_BottomBar.Pz(),weight);
139 
140  hBottomMassEta->Fill(tlv_Bottom.M(),tlv_Bottom.Eta(),weight);
141  hBottomMassEta->Fill(tlv_BottomBar.M(),tlv_BottomBar.Eta(),weight);
142 
143  hBottomMassY->Fill(tlv_Bottom.M(),tlv_Bottom.Rapidity(),weight);
144  hBottomMassY->Fill(tlv_BottomBar.M(),tlv_BottomBar.Rapidity(),weight);
145 
146  hBottomMassDeltaY->Fill(tlv_Bottom.M(),tlv_Bottom.Eta()-tlv_Bottom.Rapidity(),weight);
147  hBottomMassDeltaY->Fill(tlv_BottomBar.M(),tlv_BottomBar.Eta()-tlv_BottomBar.Rapidity(),weight);
148  }
149  if(Wplus && Wmin){
150  hWplusPz->Fill(tlv_Wplus.Pz(),weight);
151  hWminPz->Fill(tlv_Wmin.Pz(),weight);
152  }
153 }
154 
155 
156 // ------------ method called once each job just before starting event loop ------------
157 void
159 {
160  if(!dbe) return;
161  dbe->setCurrentFolder("Generator/TTbar");
162 
163  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
164 
165  hTopPt = dbe->book1D("TTbar_TopPt","t quark transverse momentum",1000,0.,1000.); hTopPt->setAxisTitle("t quark transverse momentum");
166  hTopY = dbe->book1D("TTbar_TopY","t quark rapidity",200,-5.,5.); hTopY->setAxisTitle("t quark rapidity");
167  hTopMass = dbe->book1D("TTbar_TopMass","t quark mass",500,0.,500.); hTopMass->setAxisTitle("t quark mass");
168 
169  hTTbarPt = dbe->book1D("TTbar_TTbarPt","tt pair transverse momentum",1000,0.,1000.); hTTbarPt->setAxisTitle("tt pair transverse momentum");
170  hTTbarY = dbe->book1D("TTbar_TTbarY","tt pair rapidity",200,-5.,5.); hTTbarY->setAxisTitle("tt pair rapidity");
171  hTTbarMass = dbe->book1D("TTbar_TTbarMass","tt pair mass",1000,0.,1000.); hTTbarMass->setAxisTitle("tt pair mass");
172 
173  hBottomPt = dbe->book1D("TTbar_BottomPt","b quark transverse momentum",1000,0.,1000.); hBottomPt->setAxisTitle("b quark transverse momentum");
174  hBottomEta = dbe->book1D("TTbar_BottomEta","b quark pseudo-rapidity",200,-5.,5.); hBottomEta->setAxisTitle("b quark pseudo-rapidity");
175  hBottomY = dbe->book1D("TTbar_BottomY","b quark rapidity",200,-5.,5.); hBottomY->setAxisTitle("b quark rapidity");
176  hBottomPz = dbe->book1D("TTbar_BottomPz","b quark longitudinal momentum",200,-100.,100.); hBottomPz->setAxisTitle("b quark longitudinal momentum");
177  hBottomE = dbe->book1D("TTbar_BottomE","b quark energy",1000,0.,1000.); hBottomE->setAxisTitle("b quark energy");
178  hBottomMass = dbe->book1D("TTbar_BottomMass","b quark mass",50,0.,5.); hBottomMass->setAxisTitle("b quark mass");
179 
180  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);
181  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);
182  hBottomEtaPt = dbe->book2D("TTbar_BottomEtaPt"," quark transveral momentum vs pseudorapidity",200,-5.,5.,1000,0.,1000.); hBottomEtaPt->setAxisTitle("#eta");hBottomEtaPt->setAxisTitle("P_{t} (GeV)");
183  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)");
184  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)");
185  hBottomMassEta = dbe->book2D("TTbar_BottomMassEta","b quark pseudorapidity vs mass",50,0.,5.,200,-5.,5.); hBottomMassEta->setAxisTitle("M (GeV)");hBottomMassEta->setAxisTitle("#eta");
186  hBottomMassY = dbe->book2D("TTbar_BottomMassY","b quark rapidity vs mass",50,0.,5.,200,-5.,5.); hBottomMassY->setAxisTitle("M (GeV)"); hBottomMassY->setAxisTitle("Y");
187  hBottomMassDeltaY = dbe->book2D("TTbar_BottomMassDeltaY","b quark pseudorapidity - rapidity vs mass",50,0.,50.,2000,-5.,5.); hBottomMassDeltaY->setAxisTitle("M (GeV)");hBottomMassDeltaY->setAxisTitle("Y");
188 
189  hWplusPz = dbe->book1D("TTbar_WplusPz","W+ boson longitudinal momentum",200,-100.,100.); hWplusPz->setAxisTitle("W+ boson longitudinal momentum");
190  hWminPz = dbe->book1D("TTbar_WminPz","W- boson longitudinal momentum",200,-100.,100.); hWminPz->setAxisTitle("W- boson longitudinal momentum");
191 
192 }
193 
194 // ------------ method called once each job just after ending the event loop ------------
195 void
197 {
198 }
199 
200 // ------------ method called when starting to processes a run ------------
201 void
203 {
204 }
205 
206 // ------------ method called when ending the processing of a run ------------
207 void
209 {
210 }
211 
212 // ------------ method called when starting to processes a luminosity block ------------
213 void
215 {
216 }
217 
218 // ------------ method called when ending the processing of a luminosity block ------------
219 void
221 {
222 }
223 
224 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
225 void
227  //The following says we do not know what parameters are allowed so do no validation
228  // Please change this to state exactly what you do use, even if it is no parameters
230  desc.setUnknown();
231  descriptions.addDefault(desc);
232 }
MonitorElement * hBottomYPz
MonitorElement * nEvt
TLorentzVector tlv_Bottom
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
MonitorElement * hBottomMass
MonitorElement * hBottomPtPz
virtual void endRun(edm::Run const &, edm::EventSetup const &)
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
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
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
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:845
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:429
Definition: Run.h:33
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * hBottomY