Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:11

0001 #include "Validation/EventGenerator/interface/TTbarSpinCorrHepMCAnalyzer.h"
0002 #include "Validation/EventGenerator/interface/DQMHelper.h"
0003 //
0004 // constructors and destructor
0005 //
0006 TTbarSpinCorrHepMCAnalyzer::TTbarSpinCorrHepMCAnalyzer(const edm::ParameterSet& iConfig)
0007     : genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag")),
0008       genParticlesTag_(iConfig.getParameter<edm::InputTag>("genParticlesTag")) {
0009   genEventInfoProductTagToken_ = consumes<GenEventInfoProduct>(genEventInfoProductTag_);
0010   genParticlesTagToken_ = consumes<reco::GenParticleCollection>(genParticlesTag_);
0011 }
0012 
0013 TTbarSpinCorrHepMCAnalyzer::~TTbarSpinCorrHepMCAnalyzer() {}
0014 
0015 //
0016 // member functions
0017 //
0018 
0019 // ------------ method called for each event  ------------
0020 void TTbarSpinCorrHepMCAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0021   using namespace edm;
0022 
0023   // --- the MC weights ---
0024   Handle<GenEventInfoProduct> evt_info;
0025   iEvent.getByToken(genEventInfoProductTagToken_, evt_info);
0026   if (evt_info.failedToGet())
0027     return;
0028 
0029   weight = evt_info->weight();
0030 
0031   // --- get genParticles ---
0032   Handle<reco::GenParticleCollection> genParticles;
0033   iEvent.getByToken(genParticlesTagToken_, genParticles);
0034 
0035   const reco::GenParticle* _lepton(nullptr);
0036   const reco::GenParticle* _leptonBar(nullptr);
0037 
0038   bool hasTop(false), hasTopbar(false);
0039   for (size_t i = 0; i < genParticles->size(); ++i) {
0040     const reco::GenParticle& p = (*genParticles)[i];
0041     if (p.pdgId() == 6)
0042       hasTop = true;
0043     if (p.pdgId() == -6)
0044       hasTopbar = true;
0045   }
0046 
0047   if (hasTop && hasTopbar) {
0048     // --- get status 3 leptons
0049     for (size_t i = 0; i < genParticles->size(); ++i) {
0050       const reco::GenParticle& p = (*genParticles)[i];
0051       if ((p.pdgId() == 11 || p.pdgId() == 13 || p.pdgId() == 15) && p.status() == 3) {
0052         _lepton = &p;
0053       }
0054       if ((p.pdgId() == -11 || p.pdgId() == -13 || p.pdgId() == -15) && p.status() == 3) {
0055         _leptonBar = &p;
0056       }
0057 
0058       if (_lepton && _leptonBar)
0059         break;
0060     }
0061 
0062     if (_lepton && _leptonBar) {
0063       const math::XYZTLorentzVector& lepton = _lepton->p4();
0064       const math::XYZTLorentzVector& leptonBar = _leptonBar->p4();
0065 
0066       double deltaPhi = fabs(TVector2::Phi_mpi_pi(lepton.phi() - leptonBar.phi()));
0067       _h_deltaPhi->Fill(deltaPhi, weight);
0068 
0069       double asym = (deltaPhi > CLHEP::halfpi) ? 0.5 : -0.5;
0070       _h_asym->Fill(asym, weight);
0071 
0072       math::XYZTLorentzVector llpair = lepton + leptonBar;
0073 
0074       double llpairPt = llpair.pt();
0075       _h_llpairPt->Fill(llpairPt, weight);
0076 
0077       double llpairM = llpair.M();
0078       _h_llpairM->Fill(llpairM, weight);
0079     }
0080     nEvt->Fill(0.5, weight);
0081   }
0082 }
0083 
0084 // ------------ method called once each job just before starting event loop  ------------
0085 void TTbarSpinCorrHepMCAnalyzer::bookHistograms(DQMStore::IBooker& i, edm::Run const&, edm::EventSetup const&) {
0086   ///Setting the DQM top directories
0087   std::string dir = "Generator/";
0088   dir += "TTbarSpinCorr";
0089   DQMHelper dqm(&i);
0090   i.setCurrentFolder(dir);
0091 
0092   // Number of analyzed events
0093   nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "bin", "Number of Events");
0094   _h_asym = dqm.book1dHisto("TTbar_asym", "Asymmetr", 2, -1., 1., "Asymmetry", "Number of Events");
0095   _h_deltaPhi =
0096       dqm.book1dHisto("TTbar_deltaPhi", "#Delta#phi(ll)", 320, 0, 3.2, "#Delta#phi_{ll} (rad)", "Number of Events");
0097   _h_llpairPt = dqm.book1dHisto(
0098       "TTbar_llpairPt", "Lepton pair transverse momentum", 1000, 0, 1000, "p_{T}^{ll} (GeV)", "Number of Events");
0099   _h_llpairM =
0100       dqm.book1dHisto("TTbar_llpairM", "Lepton pair invariant mass", 1000, 0, 1000, "M_{ll} (GeV)", "Number of Events");
0101 }