CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
TTbarSpinCorrHepMCAnalyzer Class Reference

#include <TTbarSpinCorrHepMCAnalyzer.h>

Inheritance diagram for TTbarSpinCorrHepMCAnalyzer:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
 TTbarSpinCorrHepMCAnalyzer (const edm::ParameterSet &)
 
 ~TTbarSpinCorrHepMCAnalyzer () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Attributes

MonitorElement_h_asym
 
MonitorElement_h_deltaPhi
 
MonitorElement_h_llpairM
 
MonitorElement_h_llpairPt
 
edm::InputTag genEventInfoProductTag_
 
edm::EDGetTokenT< GenEventInfoProductgenEventInfoProductTagToken_
 
edm::InputTag genParticlesTag_
 
edm::EDGetTokenT< reco::GenParticleCollectiongenParticlesTagToken_
 
MonitorElementnEvt
 
double weight
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 55 of file TTbarSpinCorrHepMCAnalyzer.h.

Constructor & Destructor Documentation

TTbarSpinCorrHepMCAnalyzer::TTbarSpinCorrHepMCAnalyzer ( const edm::ParameterSet iConfig)
explicit

Definition at line 6 of file TTbarSpinCorrHepMCAnalyzer.cc.

References genEventInfoProductTag_, genEventInfoProductTagToken_, genParticlesTag_, and genParticlesTagToken_.

6  :
7  genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag")),
8  genParticlesTag_(iConfig.getParameter<edm::InputTag>("genParticlesTag"))
9 {
10 
11  genEventInfoProductTagToken_=consumes<GenEventInfoProduct>(genEventInfoProductTag_);
12  genParticlesTagToken_=consumes<reco::GenParticleCollection>(genParticlesTag_);
13 
14 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesTagToken_
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoProductTagToken_
TTbarSpinCorrHepMCAnalyzer::~TTbarSpinCorrHepMCAnalyzer ( )
override

Definition at line 17 of file TTbarSpinCorrHepMCAnalyzer.cc.

17 {}

Member Function Documentation

void TTbarSpinCorrHepMCAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 24 of file TTbarSpinCorrHepMCAnalyzer.cc.

References _h_asym, _h_deltaPhi, _h_llpairM, _h_llpairPt, hiPixelPairStep_cff::deltaPhi, edm::HandleBase::failedToGet(), MonitorElement::Fill(), genEventInfoProductTagToken_, GenHFHadronMatcher_cfi::genParticles, genParticlesTagToken_, edm::Event::getByToken(), mps_fire::i, nEvt, AlCaHLTBitMon_ParallelJobs::p, reco::LeafCandidate::p4(), reco::LeafCandidate::pdgId(), reco::LeafCandidate::status(), and GenEventInfoProduct::weight().

24  {
25  using namespace edm;
26 
27  // --- the MC weights ---
29  iEvent.getByToken(genEventInfoProductTagToken_, evt_info);
30  if (evt_info.failedToGet())
31  return;
32 
33  weight = evt_info->weight() ;
34 
35  // --- get genParticles ---
37  iEvent.getByToken(genParticlesTagToken_, genParticles);
38 
39  const reco::GenParticle * _lepton (nullptr) ;
40  const reco::GenParticle * _leptonBar(nullptr) ;
41 
42  bool hasTop(false), hasTopbar(false);
43  for(size_t i = 0; i < genParticles->size(); ++ i) {
44  const reco::GenParticle & p = (*genParticles)[i];
45  if(p.pdgId() == 6) hasTop=true;
46  if(p.pdgId() == -6) hasTopbar=true;
47  }
48 
49  if(hasTop && hasTopbar){
50  // --- get status 3 leptons
51  for(size_t i = 0; i < genParticles->size(); ++ i) {
52  const reco::GenParticle & p = (*genParticles)[i];
53  if ( (p.pdgId() == 11 ||
54  p.pdgId() == 13 ||
55  p.pdgId() == 15) && p.status() == 3) { _lepton = &p ; }
56  if ( (p.pdgId() == -11 ||
57  p.pdgId() == -13 ||
58  p.pdgId() == -15) && p.status() == 3) { _leptonBar = &p ; }
59 
60  if (_lepton && _leptonBar) break;
61  }
62 
63  if (_lepton && _leptonBar) {
64 
65  const math::XYZTLorentzVector& lepton = _lepton ->p4() ;
66  const math::XYZTLorentzVector& leptonBar = _leptonBar->p4() ;
67 
68  double deltaPhi = fabs(TVector2::Phi_mpi_pi(lepton.phi() - leptonBar.phi())) ;
69  _h_deltaPhi->Fill(deltaPhi, weight) ;
70 
71  double asym = ( deltaPhi > CLHEP::halfpi ) ? 0.5 : -0.5 ;
72  _h_asym->Fill(asym, weight) ;
73 
74  math::XYZTLorentzVector llpair = lepton + leptonBar ;
75 
76  double llpairPt = llpair.pt() ;
77  _h_llpairPt->Fill(llpairPt, weight) ;
78 
79  double llpairM = llpair.M() ;
80  _h_llpairM ->Fill(llpairM , weight) ;
81 
82  }
83  nEvt->Fill(0.5 , weight) ;
84  }
85 }
int pdgId() const final
PDG identifier.
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesTagToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Definition: weight.py:1
double weight() const
void Fill(long long x)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool failedToGet() const
Definition: HandleBase.h:78
HLT enums.
int status() const final
status word
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoProductTagToken_
void TTbarSpinCorrHepMCAnalyzer::bookHistograms ( DQMStore::IBooker i,
edm::Run const &  ,
edm::EventSetup const &   
)
override

Setting the DQM top directories

Definition at line 89 of file TTbarSpinCorrHepMCAnalyzer.cc.

References _h_asym, _h_deltaPhi, _h_llpairM, _h_llpairPt, DQMHelper::book1dHisto(), dir, nEvt, DQMStore::IBooker::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

89  {
91  std::string dir="Generator/";
92  dir+="TTbarSpinCorr";
93  DQMHelper dqm(&i); i.setCurrentFolder(dir);
94 
95  // Number of analyzed events
96  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1.,"bin","Number of Events");
97  _h_asym = dqm.book1dHisto("TTbar_asym","Asymmetr", 2, -1., 1.,"Asymmetry","Number of Events");
98  _h_deltaPhi = dqm.book1dHisto("TTbar_deltaPhi","#Delta#phi(ll)", 320, 0, 3.2,"#Delta#phi_{ll} (rad)","Number of Events");
99  _h_llpairPt = dqm.book1dHisto("TTbar_llpairPt","Lepton pair transverse momentum", 1000, 0, 1000,"p_{T}^{ll} (GeV)","Number of Events");
100  _h_llpairM = dqm.book1dHisto("TTbar_llpairM","Lepton pair invariant mass", 1000, 0, 1000,"M_{ll} (GeV)","Number of Events");
101 }
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
dbl *** dir
Definition: mlp_gen.cc:35

Member Data Documentation

MonitorElement* TTbarSpinCorrHepMCAnalyzer::_h_asym
private

Definition at line 68 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* TTbarSpinCorrHepMCAnalyzer::_h_deltaPhi
private

Definition at line 69 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* TTbarSpinCorrHepMCAnalyzer::_h_llpairM
private

Definition at line 72 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* TTbarSpinCorrHepMCAnalyzer::_h_llpairPt
private

Definition at line 71 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and bookHistograms().

edm::InputTag TTbarSpinCorrHepMCAnalyzer::genEventInfoProductTag_
private

Definition at line 74 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by TTbarSpinCorrHepMCAnalyzer().

edm::EDGetTokenT<GenEventInfoProduct> TTbarSpinCorrHepMCAnalyzer::genEventInfoProductTagToken_
private

Definition at line 76 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and TTbarSpinCorrHepMCAnalyzer().

edm::InputTag TTbarSpinCorrHepMCAnalyzer::genParticlesTag_
private

Definition at line 74 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by TTbarSpinCorrHepMCAnalyzer().

edm::EDGetTokenT<reco::GenParticleCollection> TTbarSpinCorrHepMCAnalyzer::genParticlesTagToken_
private

Definition at line 77 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and TTbarSpinCorrHepMCAnalyzer().

MonitorElement* TTbarSpinCorrHepMCAnalyzer::nEvt
private

Definition at line 67 of file TTbarSpinCorrHepMCAnalyzer.h.

Referenced by analyze(), and bookHistograms().

double TTbarSpinCorrHepMCAnalyzer::weight
private

Definition at line 65 of file TTbarSpinCorrHepMCAnalyzer.h.