CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/Validation/EventGenerator/plugins/TTbar_Kinematics.cc

Go to the documentation of this file.
00001 #include "Validation/EventGenerator/interface/TTbar_Kinematics.h"
00002 #include "DQMServices/Core/interface/DQMStore.h"
00003 #include "Validation/EventGenerator/interface/PdtPdgMini.h"
00004 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00005 
00006 
00007 using namespace edm;
00008 TTbar_Kinematics::TTbar_Kinematics(const edm::ParameterSet& iConfig) :
00009   genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag"))
00010 {
00011   dbe = 0;
00012   dbe = edm::Service<DQMStore>().operator->();
00013 }
00014 
00015 
00016 TTbar_Kinematics::~TTbar_Kinematics()
00017 {
00018  
00019    // do anything here that needs to be done at desctruction time
00020    // (e.g. close files, deallocate resources etc.)
00021 
00022 }
00023 
00024 
00025 //
00026 // member functions
00027 //
00028 
00029 // ------------ method called for each event  ------------
00030 void
00031 TTbar_Kinematics::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00032 {
00033 
00034   // --- the MC weights ---
00035   Handle<GenEventInfoProduct> evt_info;
00036   iEvent.getByLabel(genEventInfoProductTag_, evt_info);
00037   if(!evt_info.isValid()) return;
00038   weight = evt_info->weight() ;
00039 
00040   // --- get TopQuarkAnalysis TtGenEvent
00041   Handle<TtGenEvent> genEvt;
00042   iEvent.getByLabel("genEvt", genEvt);
00043 
00044   if(!genEvt.isValid())return;
00045 
00046   const reco::GenParticle*        top = 0;
00047   const reco::GenParticle*    antitop = 0;
00048   const reco::GenParticle*     bottom = 0;
00049   const reco::GenParticle* antibottom = 0;
00050   const reco::GenParticle*      Wplus = 0;
00051   const reco::GenParticle*      Wmin  = 0;
00052 
00053   top        = genEvt->top();
00054   antitop    = genEvt->topBar();
00055   bottom     = genEvt->b();
00056   antibottom = genEvt->bBar();
00057   Wplus      = genEvt->wPlus();
00058   Wmin       = genEvt->wMinus();
00059 
00060   tlv_Top        = TLorentzVector(0,0,0,0) ;
00061   tlv_TopBar     = TLorentzVector(0,0,0,0) ;
00062   tlv_Bottom     = TLorentzVector(0,0,0,0) ;
00063   tlv_BottomBar  = TLorentzVector(0,0,0,0) ;
00064   tlv_Wplus      = TLorentzVector(0,0,0,0) ;
00065   tlv_Wmin       = TLorentzVector(0,0,0,0) ;
00066   tlv_TTbar      = TLorentzVector(0,0,0,0);
00067 
00068   if(top)        tlv_Top.SetPxPyPzE(top->p4().px(),top->p4().py(),top->p4().pz(),top->p4().e());
00069   if(antitop)    tlv_TopBar.SetPxPyPzE(antitop->p4().px(),antitop->p4().py(),antitop->p4().pz(),antitop->p4().e());
00070   if(bottom)     tlv_Bottom.SetPxPyPzE(bottom->p4().px(),bottom->p4().py(),bottom->p4().pz(),bottom->p4().e());
00071   if(antibottom) tlv_BottomBar.SetPxPyPzE(antibottom->p4().px(),antibottom->p4().py(),antibottom->p4().pz(),antibottom->p4().e());
00072   if(Wplus)      tlv_Wplus.SetPxPyPzE(Wplus->p4().px(),Wplus->p4().py(),Wplus->p4().pz(),Wplus->p4().e());
00073   if(Wmin)       tlv_Wmin.SetPxPyPzE(Wmin->p4().px(),Wmin->p4().py(),Wmin->p4().pz(),Wmin->p4().e());
00074   tlv_TTbar = tlv_Top + tlv_TopBar ;
00075 
00076   //---topquarkquantities---
00077   nEvt->Fill(0.5,weight);
00078   if(top && antitop){
00079     hTopPt->Fill(tlv_Top.Pt(),weight);
00080     hTopPt->Fill(tlv_TopBar.Pt(),weight);
00081     
00082     hTopY->Fill(tlv_Top.Rapidity(),weight);
00083     hTopY->Fill(tlv_TopBar.Rapidity(),weight);
00084     
00085     hTopMass->Fill(tlv_Top.M(),weight);
00086     hTopMass->Fill(tlv_TopBar.M(),weight);
00087     
00088     //---ttbarpairquantities---
00089     hTTbarPt->Fill(tlv_TTbar.Pt(),weight);
00090     hTTbarPt->Fill(tlv_TTbar.Pt(),weight);
00091     
00092     hTTbarY->Fill(tlv_TTbar.Rapidity(),weight);
00093     hTTbarY->Fill(tlv_TTbar.Rapidity(),weight);
00094     
00095     hTTbarMass->Fill(tlv_TTbar.M(),weight);
00096     hTTbarMass->Fill(tlv_TTbar.M(),weight);
00097   }
00098   if(bottom && antibottom){
00099     hBottomPt->Fill(tlv_Bottom.Pt(),weight);
00100     hBottomPt->Fill(tlv_BottomBar.Pt(),weight);
00101     
00102     hBottomEta->Fill(tlv_Bottom.Eta(),weight);
00103     hBottomEta->Fill(tlv_BottomBar.Eta(),weight);
00104     
00105     //hBottomY->Fill(math::XYZTLorentzVector(bottom->momentum()).Rapidity(),weight);
00106     //hBottomY->Fill(math::XYZTLorentzVector(antibottom->momentum()).Rapidity(),weight);
00107     
00108     hBottomY->Fill(tlv_Bottom.Rapidity(),weight);
00109     hBottomY->Fill(tlv_BottomBar.Rapidity(),weight);
00110     
00111     hBottomPz->Fill(tlv_Bottom.Pz(),weight);
00112     hBottomPz->Fill(tlv_BottomBar.Pz(),weight);
00113     
00114     hBottomE->Fill(tlv_Bottom.E(),weight);
00115     hBottomE->Fill(tlv_BottomBar.E(),weight);
00116     
00117     hBottomMass->Fill(tlv_Bottom.M(),weight);
00118     hBottomMass->Fill(tlv_BottomBar.M(),weight);
00119     
00120     hBottomPtPz->Fill(tlv_Bottom.Pt(),tlv_Bottom.Pz(),weight);
00121     hBottomPtPz->Fill(tlv_BottomBar.Pt(),tlv_BottomBar.Pz(),weight);
00122     
00123     hBottomEtaPz->Fill(tlv_Bottom.Eta(),tlv_Bottom.Pz(),weight);
00124     hBottomEtaPz->Fill(tlv_BottomBar.Eta(),tlv_BottomBar.Pz(),weight);
00125     
00126     hBottomEtaPt->Fill(tlv_Bottom.Eta(),tlv_Bottom.Pt(),weight);
00127     hBottomEtaPt->Fill(tlv_BottomBar.Eta(),tlv_BottomBar.Pt(),weight);
00128     
00129     hBottomYPz->Fill(tlv_Bottom.Rapidity(),tlv_Bottom.Pz(),weight);
00130     hBottomYPz->Fill(tlv_BottomBar.Rapidity(),tlv_BottomBar.Pz(),weight);
00131     
00132     hBottomMassPz->Fill(tlv_Bottom.M(),tlv_Bottom.Pz(),weight);
00133     hBottomMassPz->Fill(tlv_BottomBar.M(),tlv_BottomBar.Pz(),weight);
00134     
00135     hBottomMassEta->Fill(tlv_Bottom.M(),tlv_Bottom.Eta(),weight);
00136     hBottomMassEta->Fill(tlv_BottomBar.M(),tlv_BottomBar.Eta(),weight);
00137     
00138     hBottomMassY->Fill(tlv_Bottom.M(),tlv_Bottom.Rapidity(),weight);
00139     hBottomMassY->Fill(tlv_BottomBar.M(),tlv_BottomBar.Rapidity(),weight);
00140     
00141     hBottomMassDeltaY->Fill(tlv_Bottom.M(),tlv_Bottom.Eta()-tlv_Bottom.Rapidity(),weight);
00142     hBottomMassDeltaY->Fill(tlv_BottomBar.M(),tlv_BottomBar.Eta()-tlv_BottomBar.Rapidity(),weight);
00143   }
00144   if(Wplus && Wmin){
00145     hWplusPz->Fill(tlv_Wplus.Pz(),weight);
00146     hWminPz->Fill(tlv_Wmin.Pz(),weight);
00147   }
00148 }
00149 
00150 
00151 // ------------ method called once each job just before starting event loop  ------------
00152 void 
00153 TTbar_Kinematics::beginJob()
00154 {
00155   if(!dbe) return;
00156   dbe->setCurrentFolder("Generator/TTbar");
00157 
00158   nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00159 
00160   hTopPt         = dbe->book1D("TTbar_TopPt","t quark transverse momentum",1000,0.,1000.); hTopPt->setAxisTitle("t quark transverse momentum");
00161   hTopY          = dbe->book1D("TTbar_TopY","t quark rapidity",200,-5.,5.);                hTopY->setAxisTitle("t quark rapidity");
00162   hTopMass       = dbe->book1D("TTbar_TopMass","t quark mass",500,0.,500.);                hTopMass->setAxisTitle("t quark mass");
00163  
00164   hTTbarPt       = dbe->book1D("TTbar_TTbarPt","tt pair transverse momentum",1000,0.,1000.); hTTbarPt->setAxisTitle("tt pair transverse momentum");
00165   hTTbarY        = dbe->book1D("TTbar_TTbarY","tt pair rapidity",200,-5.,5.);                hTTbarY->setAxisTitle("tt pair rapidity");
00166   hTTbarMass     = dbe->book1D("TTbar_TTbarMass","tt pair mass",1000,0.,1000.);              hTTbarMass->setAxisTitle("tt pair mass");
00167 
00168   hBottomPt      = dbe->book1D("TTbar_BottomPt","b quark transverse momentum",1000,0.,1000.);     hBottomPt->setAxisTitle("b quark transverse momentum");
00169   hBottomEta     = dbe->book1D("TTbar_BottomEta","b quark pseudo-rapidity",200,-5.,5.);           hBottomEta->setAxisTitle("b quark pseudo-rapidity");
00170   hBottomY       = dbe->book1D("TTbar_BottomY","b quark rapidity",200,-5.,5.);                    hBottomY->setAxisTitle("b quark rapidity");
00171   hBottomPz      = dbe->book1D("TTbar_BottomPz","b quark longitudinal momentum",200,-100.,100.);  hBottomPz->setAxisTitle("b quark longitudinal momentum");
00172   hBottomE       = dbe->book1D("TTbar_BottomE","b quark energy",1000,0.,1000.);                   hBottomE->setAxisTitle("b quark energy");
00173   hBottomMass    = dbe->book1D("TTbar_BottomMass","b quark mass",50,0.,5.);                       hBottomMass->setAxisTitle("b quark mass");
00174   
00175   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);
00176   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);
00177   hBottomEtaPt   = dbe->book2D("TTbar_BottomEtaPt"," quark transveral   momentum vs pseudorapidity",200,-5.,5.,1000,0.,1000.);   hBottomEtaPt->setAxisTitle("#eta");hBottomEtaPt->setAxisTitle("P_{t} (GeV)");
00178   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)");
00179   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)");
00180   hBottomMassEta = dbe->book2D("TTbar_BottomMassEta","b quark pseudorapidity vs mass",50,0.,5.,200,-5.,5.);                        hBottomMassEta->setAxisTitle("M (GeV)");hBottomMassEta->setAxisTitle("#eta");
00181   hBottomMassY   = dbe->book2D("TTbar_BottomMassY","b quark rapidity vs mass",50,0.,5.,200,-5.,5.);                                hBottomMassY->setAxisTitle("M (GeV)"); hBottomMassY->setAxisTitle("Y");
00182   hBottomMassDeltaY = dbe->book2D("TTbar_BottomMassDeltaY","b quark pseudorapidity - rapidity vs mass",50,0.,50.,2000,-5.,5.);     hBottomMassDeltaY->setAxisTitle("M (GeV)");hBottomMassDeltaY->setAxisTitle("Y");
00183 
00184   hWplusPz       = dbe->book1D("TTbar_WplusPz","W+ boson longitudinal momentum",200,-100.,100.);  hWplusPz->setAxisTitle("W+ boson longitudinal momentum");
00185   hWminPz        = dbe->book1D("TTbar_WminPz","W- boson longitudinal momentum",200,-100.,100.);   hWminPz->setAxisTitle("W- boson longitudinal momentum");
00186 
00187 }
00188 
00189 // ------------ method called once each job just after ending the event loop  ------------
00190 void 
00191 TTbar_Kinematics::endJob() 
00192 {
00193 }
00194 
00195 // ------------ method called when starting to processes a run  ------------
00196 void 
00197 TTbar_Kinematics::beginRun(edm::Run const&, edm::EventSetup const&)
00198 {
00199 }
00200 
00201 // ------------ method called when ending the processing of a run  ------------
00202 void 
00203 TTbar_Kinematics::endRun(edm::Run const&, edm::EventSetup const&)
00204 {
00205 }
00206 
00207 // ------------ method called when starting to processes a luminosity block  ------------
00208 void 
00209 TTbar_Kinematics::beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
00210 {
00211 }
00212 
00213 // ------------ method called when ending the processing of a luminosity block  ------------
00214 void 
00215 TTbar_Kinematics::endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
00216 {
00217 }
00218 
00219 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
00220 void
00221 TTbar_Kinematics::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00222   //The following says we do not know what parameters are allowed so do no validation
00223   // Please change this to state exactly what you do use, even if it is no parameters
00224   edm::ParameterSetDescription desc;
00225   desc.setUnknown();
00226   descriptions.addDefault(desc);
00227 }