00001
00002
00003
00004
00005
00006
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 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 WValidation::~WValidation() {}
00030
00031 void WValidation::beginJob()
00032 {
00033 if(dbe){
00035 std::string folderName = "Generator/W";
00036 folderName+=_name;
00037 dbe->setCurrentFolder(folderName.c_str());
00038
00039
00040 nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00041
00042
00043 Wmass = dbe->book1D("Wmass","inv. Mass W", 70 ,0,140);
00044 WmassPeak = dbe->book1D("WmassPeak","inv. Mass W", 80 ,80 ,100);
00045 Wpt = dbe->book1D("Wpt","W pt",100,0,200);
00046 WptLog = dbe->book1D("WptLog","log(W pt)",100,0.,5.);
00047 Wrap = dbe->book1D("Wrap", "W y", 100, -5, 5);
00048 Wdaughters = dbe->book1D("Wdaughters", "W daughters", 60, -30, 30);
00049
00050 lepmet_mT = dbe->book1D("lepmet_mT","lepton-met transverse mass", 70 ,0,140);
00051 lepmet_mTPeak = dbe->book1D("lepmet_mTPeak","lepton-met transverse mass", 80 ,80 ,100);
00052 lepmet_pt = dbe->book1D("lepmet_pt","lepton-met",100,0,200);
00053 lepmet_ptLog = dbe->book1D("lepmet_ptLog","log(lepton-met pt)",100,0.,5.);
00054
00055 gamma_energy = dbe->book1D("gamma_energy", "photon energy in W rest frame", 200, 0., 100.);
00056 cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in W rest frame", 200, -1, 1);
00057
00058 leppt = dbe->book1D("leadpt","lepton pt", 200, 0., 200.);
00059 met = dbe->book1D("met","met", 200, 0., 200.);
00060 lepeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
00061
00062 }
00063 return;
00064 }
00065
00066 void WValidation::endJob(){return;}
00067 void WValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00068 {
00070 iSetup.getData( fPDGTable );
00071 return;
00072 }
00073 void WValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00074 void WValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00075 {
00076
00077
00078
00080 edm::Handle<HepMCProduct> evt;
00081 iEvent.getByLabel(hepmcCollection_, evt);
00082
00083
00084 const HepMC::GenEvent *myGenEvent = evt->GetEvent();
00085
00086 nEvt->Fill(0.5);
00087
00088 std::vector<const HepMC::GenParticle*> allleptons;
00089 std::vector<const HepMC::GenParticle*> allneutrinos;
00090
00091
00092 int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) ==13 ) ? 1 : 3;
00093
00094 bool vetotau = true;
00095
00096 for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
00097 if (vetotau) {
00098 if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
00099 }
00100 if((*iter)->status()==requiredstatus) {
00101
00102 if((*iter)->pdg_id()==_flavor)
00103 allleptons.push_back(*iter);
00104 else if (abs((*iter)->pdg_id()) == abs(_flavor)+1)
00105 allneutrinos.push_back(*iter);
00106 }
00107 }
00108
00109
00110 if (allleptons.empty() || allneutrinos.empty()) return;
00111
00112
00113 std::sort(allleptons.begin(), allleptons.end(), HepMCValidationHelper::sortByPt);
00114 std::sort(allneutrinos.begin(), allneutrinos.end(), HepMCValidationHelper::sortByPt);
00115
00116
00117 std::vector<const HepMC::GenParticle*> products;
00118 if (allleptons.front()->pdg_id() * allneutrinos.front()->pdg_id() > 0) return;
00119
00120
00121 if (allleptons.front()->momentum().perp() < 20. || allneutrinos.front()->momentum().perp() < 20. ) return;
00122
00123
00124 std::vector<const HepMC::GenParticle*> selectedLepton;
00125 selectedLepton.push_back(allleptons.front());
00126 std::vector<const HepMC::GenParticle*> fsrphotons;
00127 HepMCValidationHelper::findFSRPhotons(selectedLepton, myGenEvent, 0.1, fsrphotons);
00128
00129 Wdaughters->Fill(allleptons.front()->pdg_id());
00130 Wdaughters->Fill(allneutrinos.front()->pdg_id());
00131
00132
00133 TLorentzVector lep1(allleptons[0]->momentum().x(), allleptons[0]->momentum().y(), allleptons[0]->momentum().z(), allleptons[0]->momentum().t());
00134 TLorentzVector lep2(allneutrinos[0]->momentum().x(), allneutrinos[0]->momentum().y(), allneutrinos[0]->momentum().z(), allneutrinos[0]->momentum().t());
00135 TLorentzVector dilepton_mom = lep1 + lep2;
00136 TLorentzVector dilepton_andphoton_mom = dilepton_mom;
00137 std::vector<TLorentzVector> gammasMomenta;
00138 for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
00139 TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
00140 dilepton_andphoton_mom += phomom;
00141 Wdaughters->Fill(fsrphotons[ipho]->pdg_id());
00142 gammasMomenta.push_back(phomom);
00143 }
00144
00145 Wmass->Fill(dilepton_andphoton_mom.M());
00146 WmassPeak->Fill(dilepton_andphoton_mom.M());
00147 Wpt->Fill(dilepton_andphoton_mom.Pt());
00148 WptLog->Fill(log10(dilepton_andphoton_mom.Pt()));
00149 Wrap->Fill(dilepton_andphoton_mom.Rapidity());
00150
00151 TLorentzVector met_mom = HepMCValidationHelper::genMet(myGenEvent, -3., 3.);
00152 TLorentzVector lep1T(lep1.Px(), lep1.Py(), 0., lep1.Et());
00153 TLorentzVector lepmet_mom = lep1T + met_mom;
00154
00155 lepmet_mT->Fill(lepmet_mom.M());
00156 lepmet_mTPeak->Fill(lepmet_mom.M());
00157 lepmet_pt->Fill(lepmet_mom.Pt());
00158 lepmet_ptLog->Fill(log10(lepmet_mom.Pt()));
00159
00160
00161 leppt->Fill(lep1.Pt());
00162 lepeta->Fill(lep1.Eta());
00163 met->Fill(met_mom.Pt());
00164
00165
00166 TVector3 boost = dilepton_andphoton_mom.BoostVector();
00167 boost*=-1.;
00168 lep1.Boost(boost);
00169 lep2.Boost(boost);
00170 for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
00171 gammasMomenta[ipho].Boost(boost);
00172 }
00173 std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
00174
00175
00176 if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
00177 gamma_energy->Fill(gammasMomenta.front().E());
00178 double dphi = lep1.DeltaR(gammasMomenta.front());
00179 cos_theta_gamma_lepton->Fill(cos(dphi));
00180 }
00181
00182
00183 }