CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/Validation/EventGenerator/plugins/DrellYanValidation.cc

Go to the documentation of this file.
00001 /*class DrellYanValidation
00002  *  
00003  *  Class to fill dqm monitor elements from existing EDM file
00004  *
00005  *  $Date: 2011/01/24 18:27:40 $
00006  *  $Revision: 1.4 $
00007  */
00008  
00009 #include "Validation/EventGenerator/interface/DrellYanValidation.h"
00010 #include "Validation/EventGenerator/interface/HepMCValidationHelper.h"
00011 #include "TLorentzVector.h"
00012 
00013 #include "CLHEP/Units/defs.h"
00014 #include "CLHEP/Units/PhysicalConstants.h"
00015 
00016 #include "DataFormats/Math/interface/LorentzVector.h"
00017 
00018 using namespace edm;
00019 
00020 DrellYanValidation::DrellYanValidation(const edm::ParameterSet& iPSet):  
00021   hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00022   _flavor(iPSet.getParameter<int>("decaysTo")),
00023   _name(iPSet.getParameter<std::string>("name")) 
00024 {    
00025   dbe = 0;
00026   dbe = edm::Service<DQMStore>().operator->();
00027 }
00028 
00029 DrellYanValidation::~DrellYanValidation() {}
00030 
00031 void DrellYanValidation::beginJob()
00032 {
00033   if(dbe){
00035     std::string folderName = "Generator/DrellYan";
00036     folderName+=_name;
00037     dbe->setCurrentFolder(folderName.c_str());
00038     
00039     // Number of analyzed events
00040     nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00041     
00042     //Kinematics
00043     Zmass = dbe->book1D("Zmass","inv. Mass Z", 70 ,0,140);
00044     ZmassPeak = dbe->book1D("ZmassPeak","inv. Mass Z", 80 ,80 ,100);
00045     Zpt = dbe->book1D("Zpt","Z pt",100,0,200);
00046     ZptLog = dbe->book1D("ZptLog","log(Z pt)",100,0.,5.);
00047     Zrap = dbe->book1D("Zrap", "Z y", 100, -5, 5);
00048     Zdaughters = dbe->book1D("Zdaughters", "Z daughters", 60, -30, 30);
00049 
00050     dilep_mass = dbe->book1D("dilep_mass","inv. Mass dilepton", 70 ,0,140);
00051     dilep_massPeak = dbe->book1D("dilep_massPeak","inv. Mass dilepton", 80 ,80 ,100);
00052     dilep_pt = dbe->book1D("dilep_pt","dilepton pt",100,0,200);
00053     dilep_ptLog = dbe->book1D("dilep_ptLog","log(dilepton pt)",100,0.,5.);
00054     dilep_rap = dbe->book1D("dilep_rap", "dilepton y", 100, -5, 5);
00055 
00056     gamma_energy = dbe->book1D("gamma_energy", "photon energy in Z rest frame", 200, 0., 100.);
00057     cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton",      "cos_theta_gamma_lepton in Z rest frame",      200, -1, 1);
00058 
00059     leadpt = dbe->book1D("leadpt","leading lepton pt", 200, 0., 200.);    
00060     secpt  = dbe->book1D("secpt","second lepton pt", 200, 0., 200.);    
00061     leadeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
00062     seceta  = dbe->book1D("seceta","second lepton eta", 100, -5., 5.);
00063 
00064   }
00065   return;
00066 }
00067 
00068 void DrellYanValidation::endJob(){return;}
00069 void DrellYanValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00070 {
00072   iSetup.getData( fPDGTable );
00073   return;
00074 }
00075 void DrellYanValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00076 void DrellYanValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00077 { 
00078   
00079   // we *DO NOT* rely on a Z entry in the particle listings!
00080 
00082   edm::Handle<HepMCProduct> evt;
00083   iEvent.getByLabel(hepmcCollection_, evt);
00084 
00085   //Get EVENT
00086   const HepMC::GenEvent *myGenEvent = evt->GetEvent();
00087 
00088   nEvt->Fill(0.5);
00089 
00090   std::vector<const HepMC::GenParticle*> allproducts; 
00091 
00092   //requires status 1 for leptons and neutrinos (except tau)
00093   int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? 1 : 3;
00094 
00095   bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;  
00096 
00097   for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
00098     if (vetotau) {
00099       if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return; 
00100     }
00101     if((*iter)->status()==requiredstatus) {
00102       if(abs((*iter)->pdg_id())==_flavor)
00103         allproducts.push_back(*iter);
00104     }
00105   }
00106  
00107   //nothing to do if we don't have 2 particles
00108   if (allproducts.size() < 2) return; 
00109 
00110   //sort them in pt
00111   std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt); 
00112 
00113   //get the first element and the first following element with opposite charge 
00114   std::vector<const HepMC::GenParticle*> products;
00115   products.push_back(allproducts.front());
00116   const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
00117   double charge1 = PData1->charge(); 
00118   for (unsigned int i = 1; i < allproducts.size(); ++i ){
00119     const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
00120     double charge2 = PData2->charge();
00121     if (charge1*charge2 < 0) products.push_back(allproducts[i]);
00122   }
00123 
00124   //if we did not find any opposite charge pair there is nothing to do
00125   if (products.size() < 2) return; 
00126 
00127   assert(products[0]->momentum().perp() > products[1]->momentum().perp()); 
00128 
00129   //leading lepton with pt > 20.
00130   if (products[0]->momentum().perp() < 20.) return;
00131 
00132   //assemble FourMomenta
00133   TLorentzVector lep1(products[0]->momentum().x(), products[0]->momentum().y(), products[0]->momentum().z(), products[0]->momentum().t()); 
00134   TLorentzVector lep2(products[1]->momentum().x(), products[1]->momentum().y(), products[1]->momentum().z(), products[1]->momentum().t()); 
00135   TLorentzVector dilepton_mom = lep1 + lep2;
00136   TLorentzVector dilepton_andphoton_mom = dilepton_mom;
00137 
00138   //mass > 60.
00139   if (dilepton_mom.M() < 60.) return;
00140 
00141   //find possible qed fsr photons
00142   std::vector<const HepMC::GenParticle*> fsrphotons;
00143   HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
00144   
00145   Zdaughters->Fill(products[0]->pdg_id()); 
00146   Zdaughters->Fill(products[1]->pdg_id()); 
00147 
00148   std::vector<TLorentzVector> gammasMomenta;
00149   for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
00150     TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t()); 
00151     dilepton_andphoton_mom += phomom;
00152     Zdaughters->Fill(fsrphotons[ipho]->pdg_id());
00153     gammasMomenta.push_back(phomom);
00154   }  
00155   //Fill Z histograms
00156   Zmass->Fill(dilepton_andphoton_mom.M());
00157   ZmassPeak->Fill(dilepton_andphoton_mom.M());
00158   Zpt->Fill(dilepton_andphoton_mom.Pt());
00159   ZptLog->Fill(log10(dilepton_andphoton_mom.Pt())); 
00160   Zrap->Fill(dilepton_andphoton_mom.Rapidity());
00161 
00162   //Fill dilepton histograms
00163   dilep_mass->Fill(dilepton_mom.M());
00164   dilep_massPeak->Fill(dilepton_mom.M());
00165   dilep_pt->Fill(dilepton_mom.Pt());
00166   dilep_ptLog->Fill(log10(dilepton_mom.Pt()));
00167   dilep_rap->Fill(dilepton_mom.Rapidity()); 
00168 
00169   //Fill lepton histograms 
00170   leadpt->Fill(lep1.Pt());
00171   secpt->Fill(lep2.Pt());
00172   leadeta->Fill(lep1.Eta());
00173   seceta->Fill(lep2.Eta());
00174 
00175   //boost everything in the Z frame
00176   TVector3 boost = dilepton_andphoton_mom.BoostVector();
00177   boost*=-1.;
00178   lep1.Boost(boost);
00179   lep2.Boost(boost);
00180   for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
00181     gammasMomenta[ipho].Boost(boost);
00182   }
00183   std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
00184 
00185   //fill gamma histograms
00186   if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
00187     gamma_energy->Fill(gammasMomenta.front().E());
00188     double dphi = lep1.DeltaR(gammasMomenta.front()) <  lep2.DeltaR(gammasMomenta.front()) ?
00189                   lep1.DeltaPhi(gammasMomenta.front()) : lep2.DeltaPhi(gammasMomenta.front());
00190     cos_theta_gamma_lepton->Fill(cos(dphi));
00191   } 
00192 
00193 }//analyze