pythia8.cc

Go to the documentation of this file.
00001 /***********************************************************************
00002 *                                                                      *
00003 * File: pythia8_jet.cc                                                 *
00004 * This is a simple test program.                                       *
00005 * Copyright \ufffd 2005 Torbj\ufffdrn Sj\ufffdstrand                   *
00006 *                                                                      *
00007 * Modified by sj and ca                                                *
00008 *                                                                      *
00009 ***********************************************************************/
00010 #include "Pythia.h"
00011 
00012 #include "LHAFortran.h"
00013 
00014 #include "HepMCInterface.h"
00015 
00016 #include "HepMC/GenEvent.h"
00017 
00018 // this is needed to generate some random numbers                                                                                                                                                                    
00019 #include <stdlib.h>
00020 #include <time.h>
00021 
00022 #include <iostream>
00023 #include <stdio.h>
00024 #include "HepMC/IO_GenEvent.h"
00025 #include "HepMC/GenEvent.h"
00026 #include "HepMC/IO_AsciiParticles.h"
00027 #include "HepMC/SimpleVector.h"
00028 
00029 #include "TH1.h"
00030 #include "TFile.h"
00031 #include "TMath.h"
00032 #include <TCanvas.h>
00033 
00034 #include "TF1.h"
00035 #include <TH2.h>
00036 #include <TStyle.h>
00037 #include "TLegend.h"
00038 #include <TTree.h>
00039 #include <TString.h>
00040 
00041 using namespace Pythia8; 
00042 using namespace std;
00043 
00044 #include "../../include/Configuration.h"
00045 #include "../../include/JetFinder.h"
00046 #include "../../include/DiJetAnalysis.h"
00047 #include "../../include/WplusJetAnalysis.h"
00048 #include "../../include/ZAnalysis.h"
00049 #include "../../include/TauAnalysis.h"
00050 #include "../../include/TopAnalysis.h"
00051 #include "../../include/UEAnalysis.h"
00052 #include "../../include/EtMissAnalysis.h"
00053 #include "../../include/ElasScatAnalysis.h"
00054 #include "../../include/UserAnalysis.h"
00055 
00056 /**************************************************************************/
00057 
00058 int main(int argc, const char* argv[]) {
00059   
00060   if(argc <= 1 ) {
00061     std::cout << "Usage: " << argv[0] << " <Configuration File>" << std::endl;
00062     return 0;
00063   }
00064   
00065   
00066   Configuration* currentconfig = new Configuration();
00067   if(currentconfig->read(argv[1]))
00068     {
00069       std::cout << " ... skip !" << std::endl;
00070       return 0;
00071     }
00072   
00073   //argv[2] = (Pythia8DATA);
00074   
00075   JetFinder FindJet;
00076   FindJet.Init(currentconfig->max_eta,currentconfig->min_pt);
00077   FindJet.InitJetFinder(currentconfig->coneRadius,currentconfig->overlapThreshold,currentconfig->jet_ptmin,currentconfig->lepton_ptmin,currentconfig->DeltaR_lepton_track);
00078 
00079   std::vector<baseAnalysis*> analysis;
00080 
00081   // create analysis Objects and initialize them                                                                                          
00082   if (currentconfig->dijet) {
00083     std::cout << "Adding module DiJetAnalysis" << std::endl;
00084     analysis.push_back(new DiJetAnalysis);
00085   }
00086   if (currentconfig->wplusjet) {
00087     std::cout << "Adding module WplusJetAnalysis" << std::endl;
00088     analysis.push_back(new WplusJetAnalysis);
00089   }
00090   if(currentconfig->z){
00091     std::cout << "Adding module ZAnalysis" << std::endl;
00092     analysis.push_back(new ZAnalysis);
00093   }
00094   if(currentconfig->tau){
00095     std::cout << "Adding module TauAnalysis" << std::endl;
00096     analysis.push_back(new TauAnalysis);
00097   }
00098   if(currentconfig->top){
00099     std::cout << "Adding module TopAnalysis" << std::endl;
00100     analysis.push_back(new TopAnalysis);
00101   }
00102   if(currentconfig->ue){
00103     std::cout << "Adding module UEAnalysis" << std::endl;
00104     analysis.push_back(new UEAnalysis);
00105   }  
00106   if(currentconfig->etmiss){
00107     std::cout << "Adding module EtMissAnalysis" << std::endl;
00108     analysis.push_back(new EtMissAnalysis);
00109   }
00110   if(currentconfig->elasScat){
00111     std::cout << "Adding module ElasScatAnalysis" << std::endl;
00112     analysis.push_back(new ElasScatAnalysis);
00113   }
00114   if(currentconfig->user){
00115     std::cout << "Adding module UserAnalysis" << std::endl;
00116     analysis.push_back(new UserAnalysis);
00117   }
00118 
00119   // init requested analysis modules and jet finder
00120   for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00121     (*i)->Init(currentconfig->max_eta,currentconfig->min_pt);
00122     (*i)->InitJetFinder(currentconfig->coneRadius,currentconfig->overlapThreshold,currentconfig->jet_ptmin,currentconfig->lepton_ptmin,currentconfig->DeltaR_lepton_track);
00123   }
00124 
00125   
00126   HepMC::I_Pythia8 ToHepMC;
00127 
00128   // Generator. Process selection. LHC initialization.
00129   //Pythia pythia;
00130   Pythia pythia(std::getenv("Pythia8DATA"));
00131 
00132   char processline[64];
00133   // Read in commands from external file.
00134   for (unsigned int files=0; files<currentconfig->ConfigFileNames.size(); files++)
00135       {
00136           std::cout << "pythia8: Will use settings from file " << currentconfig->ConfigFileNames[files].c_str() << std::endl;
00137           pythia.readFile(currentconfig->ConfigFileNames[files].c_str());
00138       }
00139 
00140   // as done for the random seed, if it is specified
00141   pythia.readString("Random:setSeed = on");
00142   sprintf(processline,"Random:seed = %ld",currentconfig->rseed);
00143   pythia.readString(processline);
00144 
00145   // master switches for FSR, ISR and MI
00146   sprintf(processline,"PartonLevel:FSR = %d",int(currentconfig->FSR));
00147   pythia.readString(processline);
00148 
00149   sprintf(processline,"PartonLevel:ISR = %d",int(currentconfig->ISR));
00150   pythia.readString(processline);
00151 
00152   sprintf(processline,"PartonLevel:MI = %d",int(currentconfig->MI));
00153   pythia.readString(processline);
00154 
00155   // Settings can also be changed with this command
00156   //pythia.readString("WeakBosonAndParton:all = on");
00157 
00158   pythia.init( 2212, 2212, 14000.);
00159 
00160   // Shorthand for the event and for settings.
00161   Event& event = pythia.event;
00162   Settings& settings = pythia.settings;
00163 
00164   // Extract settings to be used in the main program.
00165   int nEvent = currentconfig->nevents;
00166   //sj//  int nEvent = settings.mode("Main:numberOfEvents");
00167   int nList = settings.mode("Main:numberToList");
00168   int nShow = settings.mode("Main:timesToShow");
00169   int nAbort = settings.mode("Main:timesAllowErrors");
00170   bool showChangedSettings = settings.flag("Main:showChangedSettings");
00171   bool showAllSettings = settings.flag("Main:showAllSettings");
00172   bool showAllParticleData = settings.flag("Main:showAllParticleData");
00173 
00174   // List changed data.
00175   if (showChangedSettings) settings.listChanged();
00176   if (showAllSettings) settings.listAll();
00177 
00178   // List particle data.  
00179   if (showAllParticleData) ParticleDataTable::listAll();
00180 
00181   // Begin event loop.
00182   int nShowPace = max(1,nEvent/nShow); 
00183   int iAbort = 0; 
00184   for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
00185 
00186         if (iEvent%nShowPace == 0) 
00187             cout << "pythia8: Now begin event "  << iEvent << endl;
00188 
00189     // Generate events. Quit if too many failures.
00190     if (!pythia.next()) {
00191       if (++iAbort < nAbort) continue;
00192       cout << "pythia8:  Event generation aborted prematurely, owing to error!\n"; 
00193       break;
00194     }
00195  
00196     // List first few events, both hard process and complete events.
00197     //if (iEvent < nList) { 
00198     //if (iEvent < nEvent) { 
00199     //  pythia.process.list();
00200       //event.list();
00201     // }
00202 
00203     HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
00204     ToHepMC.fill_next_event( event, hepmcevt );
00205 
00206     //hepmcevt->print();
00207 
00208     //call the Jetfinder
00209     vector<fastjet::PseudoJet> inclusive_jets = FindJet.GetJet(hepmcevt);
00210 
00211     //call analysis
00212     int ret_all = true;
00213     for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00214       (*i)->SetJet(&inclusive_jets);
00215       int ret = (*i)->Process(hepmcevt);
00216       if (ret == false) ret_all = false;
00217     }
00218     
00219     if (ret_all==false) {
00220       std::cout << "Processing event number " << iEvent << std::endl;    
00221       std::cout << "throw FAILURE"<<std::endl;
00222       //std::cout<<std::endl;
00223       //std::cout<<"HepMC Output"<<std::endl;
00224       //hepmcevt->print();
00225       std::cout<<std::endl;
00226     }
00227 
00228     //delete all the jet objects to be free for the next event
00229     FindJet.ClearEvent(hepmcevt);
00230 
00231     // clear memory from this event
00232     delete hepmcevt;
00233     // End of event loop.
00234   }
00235 
00236   // Final statistics.
00237   pythia.statistics();
00238   //  pythia.particleData.listAll();
00239   pythia.settings.writeFile("allSettings.cmnd",true);
00240 
00241   TFile* f = new TFile(currentconfig->OutputFileName.c_str(),"RECREATE");
00242   for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00243     (*i)->Finalize(f);
00244   }
00245 
00246   // clean up analysis modules
00247   for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00248     delete (*i);
00249   }
00250   analysis.resize(0);
00251 
00252   // clean up config object
00253   delete currentconfig;
00254 
00255   // Tidy up
00256   std::cout<<"pythia8: program ends"<<std::endl;
00257   f->Close();
00258 
00259   std::cout<<"Run successfully!"<<std::endl;
00260   return 0;
00261 }
00262 

Generated on Mon Jan 4 15:22:34 2010 for HepMCAnalysis by  doxygen 1.4.7