CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TTbarSpinCorrHepMCAnalyzer.cc
Go to the documentation of this file.
2 
3 //
4 // constructors and destructor
5 //
7  genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag")),
8  genParticlesTag_(iConfig.getParameter<edm::InputTag>("genParticlesTag"))
9 {
10  dbe = 0;
12 }
13 
14 
16 {
17 
18  // do anything here that needs to be done at desctruction time
19  // (e.g. close files, deallocate resources etc.)
20 
21 }
22 
23 
24 //
25 // member functions
26 //
27 
28 // ------------ method called for each event ------------
29 void
31 {
32  using namespace edm;
33 
34  // --- the MC weights ---
36  iEvent.getByLabel(genEventInfoProductTag_, evt_info);
37  if (evt_info.failedToGet())
38  return;
39 
40  weight = evt_info->weight() ;
41 
42  // --- get genParticles ---
45 
46  const reco::GenParticle * _lepton (0) ;
47  const reco::GenParticle * _leptonBar(0) ;
48 
49  bool hasTop(false), hasTopbar(false);
50  for(size_t i = 0; i < genParticles->size(); ++ i) {
51  const reco::GenParticle & p = (*genParticles)[i];
52  if(p.pdgId() == 6) hasTop=true;
53  if(p.pdgId() == -6) hasTopbar=true;
54  }
55 
56  if(hasTop && hasTopbar){
57  // --- get status 3 leptons
58  for(size_t i = 0; i < genParticles->size(); ++ i) {
59  const reco::GenParticle & p = (*genParticles)[i];
60  if ( (p.pdgId() == 11 ||
61  p.pdgId() == 13 ||
62  p.pdgId() == 15) && p.status() == 3) { _lepton = &p ; }
63  if ( (p.pdgId() == -11 ||
64  p.pdgId() == -13 ||
65  p.pdgId() == -15) && p.status() == 3) { _leptonBar = &p ; }
66 
67  if (_lepton && _leptonBar) break;
68  }
69 
70  if (_lepton && _leptonBar) {
71 
72  math::XYZTLorentzVector lepton = _lepton ->p4() ;
73  math::XYZTLorentzVector leptonBar = _leptonBar->p4() ;
74 
75  double deltaPhi = fabs(TVector2::Phi_mpi_pi(lepton.phi() - leptonBar.phi())) ;
76  _h_deltaPhi->Fill(deltaPhi, weight) ;
77 
78  double asym = ( deltaPhi > CLHEP::halfpi ) ? 0.5 : -0.5 ;
79  _h_asym->Fill(asym, weight) ;
80 
81  math::XYZTLorentzVector llpair = lepton + leptonBar ;
82 
83  double llpairPt = llpair.pt() ;
84  _h_llpairPt->Fill(llpairPt, weight) ;
85 
86  double llpairM = llpair.M() ;
87  _h_llpairM ->Fill(llpairM , weight) ;
88 
89  }
90  nEvt->Fill(0.5 , weight) ;
91  }
92 }
93 
94 
95 // ------------ method called once each job just before starting event loop ------------
96 void
98 {
99  if(dbe){
101  TString dir="Generator/";
102  dir+="TTbarSpinCorr";
103  dbe->setCurrentFolder(dir.Data());
104 
105  // Number of analyzed events
106  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
107 
108  _h_asym = dbe->book1D("TTbar_asym","Asymmetr", 2, -1., 1.);
109  _h_asym->setAxisTitle("Asymmetry");
110 
111  _h_deltaPhi = dbe->book1D("TTbar_deltaPhi","#Delta#phi(ll)", 320, 0, 3.2);
112  _h_deltaPhi->setAxisTitle("#Delta#phi(ll)");
113 
114  _h_llpairPt = dbe->book1D("TTbar_llpairPt","Lepton pair transverse momentum", 1000, 0, 1000);
115  _h_llpairPt->setAxisTitle("p_{T}(ll)");
116 
117  _h_llpairM = dbe->book1D("TTbar_llpairM","Lepton pair invariant mass", 1000, 0, 1000);
118  _h_llpairM->setAxisTitle("M(ll)");
119 
120  }
121 }
122 
123 // ------------ method called once each job just after ending the event loop ------------
124 void
126 {
127 }
128 
129 // ------------ method called when starting to processes a run ------------
130 void
132 {
133 }
134 
135 // ------------ method called when ending the processing of a run ------------
136 void
138 {
139 }
140 
141 // ------------ method called when starting to processes a luminosity block ------------
142 void
144 {
145 }
146 
147 // ------------ method called when ending the processing of a luminosity block ------------
148 void
150 {
151 }
152 
153 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
virtual void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual int status() const
status word
double Phi_mpi_pi(double x)
Definition: JetUtil.h:24
virtual void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
void Fill(long long x)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
virtual void beginRun(edm::Run const &, edm::EventSetup const &)
int iEvent
Definition: GenABIO.cc:243
virtual void analyze(const edm::Event &, const edm::EventSetup &)
TTbarSpinCorrHepMCAnalyzer(const edm::ParameterSet &)
virtual void endRun(edm::Run const &, edm::EventSetup const &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
dbl *** dir
Definition: mlp_gen.cc:35
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
Definition: Run.h:33