CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 /*class WValidation
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/WValidation.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 WValidation::WValidation(const edm::ParameterSet& iPSet): 
00021   _wmanager(iPSet),
00022   hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00023   _flavor(iPSet.getParameter<int>("decaysTo")),
00024   _name(iPSet.getParameter<std::string>("name")) 
00025 {    
00026   dbe = 0;
00027   dbe = edm::Service<DQMStore>().operator->();
00028 }
00029 
00030 WValidation::~WValidation() {}
00031 
00032 void WValidation::beginJob()
00033 {
00034   if(dbe){
00036     std::string folderName = "Generator/W";
00037     folderName+=_name;
00038     dbe->setCurrentFolder(folderName.c_str());
00039     
00040     // Number of analyzed events
00041     nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00042     
00043     //Kinematics
00044     Wmass = dbe->book1D("Wmass","inv. Mass W", 70 ,0,140);
00045     WmassPeak = dbe->book1D("WmassPeak","inv. Mass W", 80 ,80 ,100);
00046     Wpt = dbe->book1D("Wpt","W pt",100,0,200);
00047     WptLog = dbe->book1D("WptLog","log(W pt)",100,0.,5.);
00048     Wrap = dbe->book1D("Wrap", "W y", 100, -5, 5);
00049     Wdaughters = dbe->book1D("Wdaughters", "W daughters", 60, -30, 30);
00050 
00051     lepmet_mT = dbe->book1D("lepmet_mT","lepton-met transverse mass", 70 ,0,140);
00052     lepmet_mTPeak = dbe->book1D("lepmet_mTPeak","lepton-met transverse mass", 80 ,80 ,100);
00053     lepmet_pt = dbe->book1D("lepmet_pt","lepton-met",100,0,200);
00054     lepmet_ptLog = dbe->book1D("lepmet_ptLog","log(lepton-met pt)",100,0.,5.);
00055 
00056     gamma_energy = dbe->book1D("gamma_energy", "photon energy in W rest frame", 200, 0., 100.);
00057     cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton",      "cos_theta_gamma_lepton in W rest frame",      200, -1, 1);
00058 
00059     leppt = dbe->book1D("leadpt","lepton pt", 200, 0., 200.);    
00060     met   = dbe->book1D("met","met", 200, 0., 200.);    
00061     lepeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
00062 
00063   }
00064   return;
00065 }
00066 
00067 void WValidation::endJob(){return;}
00068 void WValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00069 {
00071   iSetup.getData( fPDGTable );
00072   return;
00073 }
00074 void WValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00075 void WValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00076 { 
00077   
00078   // we *DO NOT* rely on a Z entry in the particle listings!
00079 
00081   edm::Handle<HepMCProduct> evt;
00082   iEvent.getByLabel(hepmcCollection_, evt);
00083 
00084   //Get EVENT
00085   const HepMC::GenEvent *myGenEvent = evt->GetEvent();
00086 
00087   double weight = _wmanager.weight(iEvent);
00088 
00089   nEvt->Fill(0.5,weight);
00090 
00091   std::vector<const HepMC::GenParticle*> allleptons; 
00092   std::vector<const HepMC::GenParticle*> allneutrinos; 
00093 
00094   //requires status 1 for leptons and neutrinos (except tau)
00095   int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) ==13 ) ? 1 : 3;
00096 
00097   bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;  
00098 
00099   for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
00100     if (vetotau) {
00101       if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
00102     }
00103     if((*iter)->status()==requiredstatus) {
00104       //@todo: improve this selection   
00105       if((*iter)->pdg_id()==_flavor)
00106         allleptons.push_back(*iter);
00107       else if (abs((*iter)->pdg_id()) == abs(_flavor)+1)
00108         allneutrinos.push_back(*iter);  
00109     }
00110   }
00111  
00112   //nothing to do if we don't have 2 particles
00113   if (allleptons.empty() || allneutrinos.empty()) return;
00114 
00115   //sort them in pt
00116   std::sort(allleptons.begin(), allleptons.end(), HepMCValidationHelper::sortByPt); 
00117   std::sort(allneutrinos.begin(), allneutrinos.end(), HepMCValidationHelper::sortByPt); 
00118 
00119   //get the first lepton and the first neutrino, and check that one is particle one is antiparticle (product of pdgids < 0) 
00120   std::vector<const HepMC::GenParticle*> products;
00121   if (allleptons.front()->pdg_id() * allneutrinos.front()->pdg_id() > 0) return;        
00122 
00123   //require at least 20 GeV on the lepton
00124   if (allleptons.front()->momentum().perp() < 20. || allneutrinos.front()->momentum().perp() < 20. ) return;
00125 
00126   //find possible qed fsr photons
00127   std::vector<const HepMC::GenParticle*> selectedLepton;
00128   selectedLepton.push_back(allleptons.front()); 
00129   std::vector<const HepMC::GenParticle*> fsrphotons;
00130   HepMCValidationHelper::findFSRPhotons(selectedLepton, myGenEvent, 0.1, fsrphotons);
00131 
00132   Wdaughters->Fill(allleptons.front()->pdg_id(),weight); 
00133   Wdaughters->Fill(allneutrinos.front()->pdg_id(),weight); 
00134  
00135   //assemble FourMomenta
00136   TLorentzVector lep1(allleptons[0]->momentum().x(), allleptons[0]->momentum().y(), allleptons[0]->momentum().z(), allleptons[0]->momentum().t()); 
00137   TLorentzVector lep2(allneutrinos[0]->momentum().x(), allneutrinos[0]->momentum().y(), allneutrinos[0]->momentum().z(), allneutrinos[0]->momentum().t()); 
00138   TLorentzVector dilepton_mom = lep1 + lep2;
00139   TLorentzVector dilepton_andphoton_mom = dilepton_mom;
00140   std::vector<TLorentzVector> gammasMomenta;
00141   for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
00142     TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t()); 
00143     dilepton_andphoton_mom += phomom;
00144     Wdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
00145     gammasMomenta.push_back(phomom);
00146   }  
00147   //Fill "true" W histograms
00148   Wmass->Fill(dilepton_andphoton_mom.M(),weight);
00149   WmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
00150   Wpt->Fill(dilepton_andphoton_mom.Pt(),weight);
00151   WptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight); 
00152   Wrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
00153 
00154   TLorentzVector met_mom = HepMCValidationHelper::genMet(myGenEvent, -3., 3.);
00155   TLorentzVector lep1T(lep1.Px(), lep1.Py(), 0., lep1.Et());
00156   TLorentzVector lepmet_mom = lep1T + met_mom;
00157   //Fill lepmet histograms
00158   lepmet_mT->Fill(lepmet_mom.M(),weight);
00159   lepmet_mTPeak->Fill(lepmet_mom.M(),weight);
00160   lepmet_pt->Fill(lepmet_mom.Pt(),weight);
00161   lepmet_ptLog->Fill(log10(lepmet_mom.Pt()),weight);
00162 
00163   //Fill lepton histograms 
00164   leppt->Fill(lep1.Pt(),weight);
00165   lepeta->Fill(lep1.Eta(),weight);
00166   met->Fill(met_mom.Pt(),weight);       
00167 
00168   //boost everything in the W frame
00169   TVector3 boost = dilepton_andphoton_mom.BoostVector();
00170   boost*=-1.;
00171   lep1.Boost(boost);
00172   lep2.Boost(boost);
00173   for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
00174     gammasMomenta[ipho].Boost(boost);
00175   }
00176   std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
00177 
00178   //fill gamma histograms
00179   if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
00180     gamma_energy->Fill(gammasMomenta.front().E(),weight);
00181     double dphi = lep1.DeltaR(gammasMomenta.front());
00182     cos_theta_gamma_lepton->Fill(cos(dphi),weight);
00183   } 
00184 
00185 
00186 }//analyze