HepMC Reference Documentation

HepMC

main31.cc

Go to the documentation of this file.
00001 // main31.cc is a part of the PYTHIA event generator.
00002 // Copyright (C) 2011 Mikhail Kirsanov, Torbjorn Sjostrand.
00003 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
00004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
00005 
00006 // This is a simple test program.
00007 // It illustrates how HepMC can be interfaced to Pythia8.
00008 // It studies the charged multiplicity distribution at the LHC.
00009 // HepMC events are output to the hepmcout31.dat file.
00010 // Written by Mikhail Kirsanov based on main01.cc.
00011 
00012 // WARNING: typically one needs 25 MB/100 events at the LHC.
00013 // Therefore large event samples may be impractical.
00014 
00015 #include "Pythia.h"
00016 #include "HepMCInterface.h"
00017 
00018 #include "HepMC/GenEvent.h"
00019 #include "HepMC/IO_GenEvent.h"
00020 
00021 // Following line is a deprecated alternative, removed in recent versions
00022 //#include "HepMC/IO_Ascii.h"
00023 //#include "HepMC/IO_AsciiParticles.h"
00024 
00025 // Following line to be used with HepMC 2.04 onwards.
00026 #ifdef HEPMC_HAS_UNITS
00027 #include "HepMC/Units.h"
00028 #endif
00029 
00030 using namespace Pythia8; 
00031 
00032 int main() {
00033 
00034   // Interface for conversion from Pythia8::Event to HepMC one. 
00035   HepMC::I_Pythia8 ToHepMC;
00036   //  ToHepMC.set_crash_on_problem();
00037 
00038   // Specify file where HepMC events will be stored.
00039   HepMC::IO_GenEvent ascii_io("hepmcout31.dat", std::ios::out);
00040   // Following line is a deprecated alternative, removed in recent versions
00041   // HepMC::IO_Ascii ascii_io("hepmcout31.dat", std::ios::out);
00042   // Line below is an eye-readable one-way output, uncomment the include above
00043   // HepMC::IO_AsciiParticles ascii_io("hepmcout31.dat", std::ios::out);
00044 
00045   // Generator. Process selection. LHC initialization. Histogram.
00046   Pythia pythia;
00047   pythia.readString("HardQCD:all = on");    
00048   pythia.readString("PhaseSpace:pTHatMin = 20.");    
00049   pythia.init( 2212, 2212, 14000.);
00050   Hist mult("charged multiplicity", 100, -0.5, 799.5);
00051 
00052   // Begin event loop. Generate event. Skip if error. List first one.
00053   for (int iEvent = 0; iEvent < 100; ++iEvent) {
00054     if (!pythia.next()) continue;
00055     if (iEvent < 1) {pythia.info.list(); pythia.event.list();} 
00056 
00057     // Find number of all final charged particles and fill histogram.
00058     int nCharged = 0;
00059     for (int i = 0; i < pythia.event.size(); ++i) 
00060       if (pythia.event[i].isFinal() && pythia.event[i].isCharged()) 
00061         ++nCharged; 
00062     mult.fill( nCharged );
00063 
00064     // Construct new empty HepMC event. 
00065 #ifdef HEPMC_HAS_UNITS
00066     // This form with arguments is only meaningful for HepMC 2.04 onwards, 
00067     // and even then unnecessary if HepMC was built with GeV and mm as units..
00068     HepMC::GenEvent* hepmcevt = new HepMC::GenEvent( 
00069       HepMC::Units::GEV, HepMC::Units::MM);
00070 #else
00071     // This form is needed for backwards compatibility. 
00072     // In HepMCInterface.cc a conversion from GeV to MeV will be done.
00073     HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
00074 #endif
00075 
00076     // Fill HepMC event, including PDF info.
00077     ToHepMC.fill_next_event( pythia, hepmcevt );
00078     // This alternative older method fills event, without PDF info.
00079     // ToHepMC.fill_next_event( pythia.event, hepmcevt );
00080 
00081     // Write the HepMC event to file. Done with it.
00082     ascii_io << hepmcevt;
00083     delete hepmcevt;
00084 
00085   // End of event loop. Statistics. Histogram. 
00086   }
00087   pythia.statistics();
00088   cout << mult; 
00089 
00090   // Done.
00091   return 0;
00092 }

Generated on Fri Feb 17 00:31:25 2012 for HepMC by  doxygen 1.4.7