00001
00002
00003
00004
00005
00006
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 _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 DrellYanValidation::~DrellYanValidation() {}
00031
00032 void DrellYanValidation::beginJob()
00033 {
00034 if(dbe){
00036 std::string folderName = "Generator/DrellYan";
00037 folderName+=_name;
00038 dbe->setCurrentFolder(folderName.c_str());
00039
00040
00041 nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00042
00043
00044 Zmass = dbe->book1D("Zmass","inv. Mass Z", 70 ,0,140);
00045 ZmassPeak = dbe->book1D("ZmassPeak","inv. Mass Z", 80 ,80 ,100);
00046 Zpt = dbe->book1D("Zpt","Z pt",100,0,200);
00047 ZptLog = dbe->book1D("ZptLog","log(Z pt)",100,0.,5.);
00048 Zrap = dbe->book1D("Zrap", "Z y", 100, -5, 5);
00049 Zdaughters = dbe->book1D("Zdaughters", "Z daughters", 60, -30, 30);
00050
00051 dilep_mass = dbe->book1D("dilep_mass","inv. Mass dilepton", 70 ,0,140);
00052 dilep_massPeak = dbe->book1D("dilep_massPeak","inv. Mass dilepton", 80 ,80 ,100);
00053 dilep_pt = dbe->book1D("dilep_pt","dilepton pt",100,0,200);
00054 dilep_ptLog = dbe->book1D("dilep_ptLog","log(dilepton pt)",100,0.,5.);
00055 dilep_rap = dbe->book1D("dilep_rap", "dilepton y", 100, -5, 5);
00056
00057 gamma_energy = dbe->book1D("gamma_energy", "photon energy in Z rest frame", 200, 0., 100.);
00058 cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in Z rest frame", 200, -1, 1);
00059
00060 leadpt = dbe->book1D("leadpt","leading lepton pt", 200, 0., 200.);
00061 secpt = dbe->book1D("secpt","second lepton pt", 200, 0., 200.);
00062 leadeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
00063 seceta = dbe->book1D("seceta","second lepton eta", 100, -5., 5.);
00064
00065 }
00066 return;
00067 }
00068
00069 void DrellYanValidation::endJob(){return;}
00070 void DrellYanValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00071 {
00073 iSetup.getData( fPDGTable );
00074 return;
00075 }
00076 void DrellYanValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00077 void DrellYanValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00078 {
00079
00080
00081
00083 edm::Handle<HepMCProduct> evt;
00084 iEvent.getByLabel(hepmcCollection_, evt);
00085
00086
00087 const HepMC::GenEvent *myGenEvent = evt->GetEvent();
00088
00089 double weight = _wmanager.weight(iEvent);
00090
00091
00092
00093 nEvt->Fill(0.5,weight);
00094
00095 std::vector<const HepMC::GenParticle*> allproducts;
00096
00097
00098 int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? 1 : 3;
00099
00100 bool vetotau = true;
00101
00102 for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
00103 if (vetotau) {
00104 if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
00105 }
00106 if((*iter)->status()==requiredstatus) {
00107 if(abs((*iter)->pdg_id())==_flavor)
00108 allproducts.push_back(*iter);
00109 }
00110 }
00111
00112
00113 if (allproducts.size() < 2) return;
00114
00115
00116 std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt);
00117
00118
00119 std::vector<const HepMC::GenParticle*> products;
00120 products.push_back(allproducts.front());
00121 const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
00122 double charge1 = PData1->charge();
00123 for (unsigned int i = 1; i < allproducts.size(); ++i ){
00124 const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
00125 double charge2 = PData2->charge();
00126 if (charge1*charge2 < 0) products.push_back(allproducts[i]);
00127 }
00128
00129
00130 if (products.size() < 2) return;
00131
00132 assert(products[0]->momentum().perp() >= products[1]->momentum().perp());
00133
00134
00135 if (products[0]->momentum().perp() < 20.) return;
00136
00137
00138 TLorentzVector lep1(products[0]->momentum().x(), products[0]->momentum().y(), products[0]->momentum().z(), products[0]->momentum().t());
00139 TLorentzVector lep2(products[1]->momentum().x(), products[1]->momentum().y(), products[1]->momentum().z(), products[1]->momentum().t());
00140 TLorentzVector dilepton_mom = lep1 + lep2;
00141 TLorentzVector dilepton_andphoton_mom = dilepton_mom;
00142
00143
00144 if (dilepton_mom.M() < 60.) return;
00145
00146
00147 std::vector<const HepMC::GenParticle*> fsrphotons;
00148 HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
00149
00150 Zdaughters->Fill(products[0]->pdg_id(),weight);
00151 Zdaughters->Fill(products[1]->pdg_id(),weight);
00152
00153 std::vector<TLorentzVector> gammasMomenta;
00154 for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
00155 TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
00156 dilepton_andphoton_mom += phomom;
00157 Zdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
00158 gammasMomenta.push_back(phomom);
00159 }
00160
00161 Zmass->Fill(dilepton_andphoton_mom.M(),weight);
00162 ZmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
00163 Zpt->Fill(dilepton_andphoton_mom.Pt(),weight);
00164 ZptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight);
00165 Zrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
00166
00167
00168 dilep_mass->Fill(dilepton_mom.M(),weight);
00169 dilep_massPeak->Fill(dilepton_mom.M(),weight);
00170 dilep_pt->Fill(dilepton_mom.Pt(),weight);
00171 dilep_ptLog->Fill(log10(dilepton_mom.Pt()),weight);
00172 dilep_rap->Fill(dilepton_mom.Rapidity(),weight);
00173
00174
00175 leadpt->Fill(lep1.Pt(),weight);
00176 secpt->Fill(lep2.Pt(),weight);
00177 leadeta->Fill(lep1.Eta(),weight);
00178 seceta->Fill(lep2.Eta(),weight);
00179
00180
00181 TVector3 boost = dilepton_andphoton_mom.BoostVector();
00182 boost*=-1.;
00183 lep1.Boost(boost);
00184 lep2.Boost(boost);
00185 for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
00186 gammasMomenta[ipho].Boost(boost);
00187 }
00188 std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
00189
00190
00191 if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
00192 gamma_energy->Fill(gammasMomenta.front().E(),weight);
00193 double dphi = lep1.DeltaR(gammasMomenta.front()) < lep2.DeltaR(gammasMomenta.front()) ?
00194 lep1.DeltaPhi(gammasMomenta.front()) : lep2.DeltaPhi(gammasMomenta.front());
00195 cos_theta_gamma_lepton->Fill(cos(dphi),weight);
00196 }
00197
00198 }