00001
00002
00003 #include "ThePEG/Repository/EventGenerator.h"
00004 #include "ThePEG/Repository/Repository.h"
00005 #include "ThePEG/Persistency/PersistentIStream.h"
00006 #include "ThePEG/Utilities/DynamicLoader.h"
00007 #include "ThePEG/EventRecord/Event.h"
00008 #include "ThePEG/EventRecord/SubProcess.h"
00009 #include "ThePEG/Handlers/XComb.h"
00010 #include "ThePEG/Handlers/EventHandler.h"
00011 #include "ThePEG/PDF/PartonExtractor.h"
00012 #include "ThePEG/PDF/PDF.h"
00013 #include "ThePEG/Vectors/HepMCConverter.h"
00014 #include "ThePEGHepMC_1.3.0.h"
00015
00016 #include "HepMC/GenEvent.h"
00017 #include "HepMC/IO_GenEvent.h"
00018 #include "HepMC/GenEvent.h"
00019 #include "HepMC/IO_AsciiParticles.h"
00020 #include "HepMC/SimpleVector.h"
00021
00022 #include <iostream>
00023 #include <stdio.h>
00024 #include <string>
00025
00026 #include "TH1.h"
00027 #include "TFile.h"
00028 #include "TMath.h"
00029 #include <TCanvas.h>
00030
00031 #include "TF1.h"
00032 #include <TH2.h>
00033 #include <TStyle.h>
00034 #include "TLegend.h"
00035 #include <TTree.h>
00036 #include <TString.h>
00037
00038 using namespace ThePEG;
00039
00040
00041 #include "../../include/Configuration.h"
00042 #include "../../include/JetFinder.h"
00043 #include "../../include/TopAnalysis.h"
00044 #include "../../include/TauAnalysis.h"
00045 #include "../../include/DiJetAnalysis.h"
00046 #include "../../include/WplusJetAnalysis.h"
00047 #include "../../include/ZAnalysis.h"
00048 #include "../../include/UEAnalysis.h"
00049 #include "../../include/EtMissAnalysis.h"
00050 #include "../../include/ElasScatAnalysis.h"
00051 #include "../../include/UserAnalysis.h"
00052 #include "readHerwigConfigFile.h"
00053
00054 int main(int argc, const char* argv[]) {
00055
00056 if(argc <= 1 ) {
00057 std::cout << "Usage: " << argv[0] << " <Configuration File>" << std::endl;
00058 return 0;
00059 }
00060
00061
00062 Configuration* currentconfig = new Configuration();
00063 if(currentconfig->read(argv[1]))
00064 {
00065 std::cout << " ... skip !" << std::endl;
00066 return 0;
00067 }
00068
00069 JetFinder FindJet;
00070 FindJet.Init(currentconfig->max_eta,currentconfig->min_pt);
00071 FindJet.InitJetFinder(currentconfig->coneRadius,currentconfig->overlapThreshold,currentconfig->jet_ptmin,currentconfig->lepton_ptmin,currentconfig->DeltaR_lepton_track);
00072
00073
00074 std::vector<baseAnalysis*> analysis;
00075
00076
00077 if (currentconfig->dijet) {
00078 std::cout << "Adding module DiJetAnalysis" << std::endl;
00079 analysis.push_back(new DiJetAnalysis);
00080 }
00081 if (currentconfig->wplusjet) {
00082 std::cout << "Adding module WplusJetAnalysis" << std::endl;
00083 analysis.push_back(new WplusJetAnalysis);
00084 }
00085 if(currentconfig->z){
00086 std::cout << "Adding module ZAnalysis" << std::endl;
00087 analysis.push_back(new ZAnalysis);
00088 }
00089 if(currentconfig->tau){
00090 std::cout << "Adding module TauAnalysis" << std::endl;
00091 analysis.push_back(new TauAnalysis);
00092 }
00093 if(currentconfig->top){
00094 std::cout << "Adding module TopAnalysis" << std::endl;
00095 analysis.push_back(new TopAnalysis);
00096 }
00097 if(currentconfig->ue){
00098 std::cout << "Adding module UEAnalysis" << std::endl;
00099 analysis.push_back(new UEAnalysis);
00100 }
00101 if(currentconfig->etmiss){
00102 std::cout << "Adding module EtMissAnalysis" << std::endl;
00103 analysis.push_back(new EtMissAnalysis);
00104 }
00105 if(currentconfig->elasScat){
00106 std::cout << "Adding module ElasScatAnalysis" << std::endl;
00107 analysis.push_back(new ElasScatAnalysis);
00108 }
00109 if(currentconfig->user){
00110 std::cout << "Adding module UserAnalysis" << std::endl;
00111 analysis.push_back(new UserAnalysis);
00112 }
00113
00114
00115
00116 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00117 (*i)->Init(currentconfig->max_eta,currentconfig->min_pt);
00118 (*i)->InitJetFinder(currentconfig->coneRadius,currentconfig->overlapThreshold,currentconfig->jet_ptmin,currentconfig->lepton_ptmin,currentconfig->DeltaR_lepton_track);
00119 }
00120
00121
00122
00123 std::cout<<"Add a path to seach for dynamically linkable libraries..."<<std::endl;
00124 ThePEG::DynamicLoader::appendPath(HWMODULEDIR);
00125 ThePEG::DynamicLoader::appendPath(THEPEGMODULEDIR);
00126
00127
00128 std::cout<<"Load repository from a given file..."<<std::endl;
00129 const string repopath = string(HWREPODIR) + "/HerwigDefaults.rpo";
00130 ThePEG::Repository::load(repopath);
00131
00132
00133 std::cout<<"Setup process,interpret the command in cmd and return possible messages..."<<std::endl;
00134
00135
00136 char exec_command[100];
00137 sprintf(exec_command,"set LHCGenerator:RandomNumberGenerator:Seed %ld",currentconfig->rseed);
00138 ThePEG::Repository::exec("cd /Herwig/Generators", cout);
00139 ThePEG::Repository::exec(exec_command, cout);
00140
00141
00142 std::vector< std::string > herwig_commands;
00143
00144 for (unsigned int files=0; files<(currentconfig->ConfigFileNames).size(); files++)
00145 readConfigFile(&herwig_commands,(currentconfig->ConfigFileNames)[files]);
00146
00147
00148 for(unsigned int i=0; i<herwig_commands.size(); i++)
00149 ThePEG::Repository::exec( herwig_commands[i].c_str(), cout);
00150
00151
00152 ThePEG::Repository::exec("cd /Herwig/Shower", cout);
00153
00154 sprintf(exec_command,"set SplittingGenerator:ISR %c",currentconfig->ISR);
00155 ThePEG::Repository::exec(exec_command, cout);
00156
00157 sprintf(exec_command,"set SplittingGenerator:FSR %c",currentconfig->FSR);
00158 ThePEG::Repository::exec(exec_command, cout);
00159
00160 sprintf(exec_command,"set ShowerHandler:MPI %c",currentconfig->MI);
00161 ThePEG::Repository::exec(exec_command, cout);
00162
00163
00164
00165 ThePEG::Repository::update();
00166
00167
00168 ThePEG::EGPtr m_hw;
00169
00170
00171 ThePEG::EventPtr m_event;
00172
00173
00174 std::cout<<"make a run object ..."<<std::endl;
00175 ThePEG::EGPtr tmpEG = ThePEG::BaseRepository::GetObject<ThePEG::EGPtr>("/Herwig/Generators/LHCGenerator");
00176 m_hw = ThePEG::Repository::makeRun(tmpEG, "LHC");
00177
00178
00179 std::cout<<"Herwig++ initializing..."<<std::endl;
00180 m_hw->initialize();
00181
00182
00183 std::cout<<"Generating events and filling histograms..."<<std::endl;
00184
00185
00186 for (unsigned int iEvent = 0; iEvent < currentconfig->nevents; ++iEvent) {
00187
00188
00189 assert(m_hw);
00190 m_event = m_hw->shoot();
00191
00192
00193
00194 HepMC::GenEvent* evt = new HepMC::GenEvent();
00195 ThePEG::HepMCConverter<HepMC::GenEvent>::convert(*m_event, *evt, false, ThePEG::GeV);
00196
00197
00198 std::cout << "Processing event number " << iEvent << std::endl;
00199 if (iEvent < 0) evt->print();
00200
00201
00202 vector<fastjet::PseudoJet> inclusive_jets = FindJet.GetJet(evt);
00203
00204
00205 int ret_all = true;
00206 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00207 (*i)->SetJet(&inclusive_jets);
00208 int ret = (*i)->Process(evt);
00209 if (ret == false) ret_all = false;
00210 }
00211
00212 if (ret_all==false) {
00213 std::cout << "Processing event number " << iEvent << std::endl;
00214 std::cout << "throw FAILURE"<<std::endl;
00215
00216
00217 std::cout<<std::endl;
00218 }
00219
00220
00221 FindJet.ClearEvent(evt);
00222
00223
00224 delete evt;
00225 }
00226
00227 TFile* f = new TFile(currentconfig->OutputFileName.c_str(),"RECREATE");
00228 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00229 (*i)->Finalize(f);
00230 }
00231
00232
00233 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00234 delete (*i);
00235 }
00236 analysis.resize(0);
00237
00238
00239 delete currentconfig;
00240
00241
00242 std::cout<<"Herwig++ finalizing..."<<std::endl;
00243 assert(m_hw);
00244 m_hw->finalize();
00245 ThePEG::Repository::cleanup();
00246
00247 std::cout<<"Run successfully!"<<std::endl;
00248 return 0;
00249 }