hepmc

Subversion Repositories:
Compare Path: Rev
With Path: Rev
/ @ 508  →  / @ 509
/trunk/ChangeLog
@@ -2,6 +2,8 @@
2012-01-05 Lynn Garren
 
* src/GenEvent.cc: fix error text for GenEvent::use_length_unit
* examples/pythia8: add main31 and main32 from pythia8
 
2012-01-05 Lynn Garren
 
/trunk/configure.ac
@@ -163,7 +163,10 @@
test/testMultipleCopies.cc
test/testStreamIO.cc
examples/GNUmakefile.example
examples/fio/GNUmakefile.example])
examples/fio/GNUmakefile.example
examples/pythia8/config.csh
examples/pythia8/config.sh
examples/pythia8/GNUmakefile.example])
 
AC_CONFIG_FILES([test/testHepMC.sh], [chmod +x test/testHepMC.sh])
AC_CONFIG_FILES([test/testFlow.sh], [chmod +x test/testFlow.sh])
/trunk/Makefile.am
@@ -2,6 +2,6 @@
 
includedir = $(prefix)/include
 
SUBDIRS = HepMC src fio test examples examples/fio doc
SUBDIRS = HepMC src fio test examples examples/fio examples/pythia8 doc
# list all subdirectories - for distribution and cleaning
DIST_SUBDIRS = HepMC src fio test examples examples/fio doc
DIST_SUBDIRS = HepMC src fio test examples examples/fio examples/pythia8 doc
New file
/trunk/examples/pythia8/config.sh.in
@@ -0,0 +1,8 @@
#!/bin/sh
if [ ! $?LD_LIBRARY_PATH ]; then
export LD_LIBRARY_PATH=@prefix@/lib
fi
if [ $?LD_LIBRARY_PATH ]; then
export LD_LIBRARY_PATH=${LD_LIBRARY_PATH}:@prefix@/lib
fi
export PYTHIA8DATA=${PYTHIA8_HOME}/xmldoc
New file
/trunk/examples/pythia8/GNUmakefile.example.in
@@ -0,0 +1,74 @@
################################################################################
# Makefile for building a HepMC PYTHIA8 example
#
# Lynn Garren garren@fnal.gov
#
# main31 and main32 have been copied from PYTHIA 8.157
# you must source either config.csh or config.sh first
#
# ./main31.exe
# ./main32.exe main32.cmnd hepmcout32.dat
#
################################################################################ Define directory paths
# You may have to change PYTHIA8_HOME and/or other variables
#
HepMCdir = @prefix@
HEPMCLOCATION = $(HepMCdir)
HepMClib = -L$(HEPMCLOCATION)/lib -Wl,-rpath -Wl,$(HEPMCLOCATION)/lib \
-lHepMC
Pythia_LIB = -L$(PYTHIA8_HOME)/lib -Wl,-rpath -Wl,$(PYTHIA8_HOME)/lib \
-lpythia8 -llhapdfdummy -lhepmcinterface
 
################################################################################ Compiler options
#
CXX = @CXX@
INCLUDES = -I$(HEPMCLOCATION)/include -I$(PYTHIA8_HOME)/include
CXXFLAGS = @AM_CXXFLAGS@ @CXXFLAGS@ -Wshadow -fbounds-check $(INCLUDES)
ifeq "$(CXX)" "g++"
FLAGS = $(DFLG) -fno-second-underscore $(INCDIR)
else
FLAGS = $(DFLG) $(INCDIR)
endif
 
LINK_LIBS = @LDFLAGS@
 
UNAME = $(shell uname)
ifneq "$(UNAME)" "Darwin"
LINK_LIBS += -Wl,-rpath -Wl,$(HepMCdir)/lib
endif
 
HDRS = $(HepMCdir)/include/HepMC/*.h *.h
EXAMPLES = main31.exe main32.exe
 
###############################################################################
#
.SUFFIXES: .o .cc .exe
 
all: $(EXAMPLES)
 
main31.exe: main31.o
$(CXX) $(CXXFLAGS) -o $@ main31.o \
$(Pythia_LIB) \
$(HepMClib)
 
main32.exe: main32.o
$(CXX) $(CXXFLAGS) -o $@ main32.o \
$(Pythia_LIB) \
$(HepMClib)
 
###############################################################################
#
.cc.o: $(HDRS) $<
$(CXX) $(CXXFLAGS) -c $< -o $@
 
###############################################################################
# gmake clean removes all garbage from HepMC directories.
# gmake distclean removes all compiled libraries, executables, +garbage
#
clean:
rm -f *.o
 
distclean:
$(MAKE) clean --no-print-directory
rm -f $(EXAMPLES)
rm -f *.dat
New file
/trunk/examples/pythia8/config.csh.in
@@ -0,0 +1,7 @@
#!/bin/csh
if( ! $?LD_LIBRARY_PATH ) then
setenv LD_LIBRARY_PATH @prefix@/lib
else
setenv LD_LIBRARY_PATH ${LD_LIBRARY_PATH}:@prefix@/lib
endif
setenv PYTHIA8DATA ${PYTHIA8_HOME}/xmldoc
New file
/trunk/examples/pythia8/Makefile.am
@@ -0,0 +1,18 @@
# This makefile is used to copy files during the "make install" step
 
INSTALLDIR = $(pkgdatadir)/examples/pythia8
 
# files to distribute
EXTRA_DIST = \
main31.cc \
main32.cc \
main32.cmnd \
README
 
install-data-local:
$(mkinstalldirs) $(DESTDIR)$(INSTALLDIR)
$(INSTALL_DATA) GNUmakefile.example $(DESTDIR)$(INSTALLDIR)/GNUmakefile
$(INSTALL_DATA) config.csh $(DESTDIR)$(INSTALLDIR)/config.csh
$(INSTALL_DATA) config.sh $(DESTDIR)$(INSTALLDIR)/config.sh
for file in $(EXTRA_DIST); do \
$(INSTALL_DATA) $(srcdir)/$$file $(DESTDIR)$(INSTALLDIR)/$$file; done
New file
/trunk/examples/pythia8/main32.cmnd
@@ -0,0 +1,32 @@
! File: main32.cmnd
! This file contains commands to be read in for a Pythia8 run.
! Lines not beginning with a letter or digit are comments.
! Names are case-insensitive - but spellings-sensitive!
! The changes here are illustrative, not always physics-motivated.
 
! 1) Settings that will be used in a main program.
Main:numberOfEvents = 200 ! number of events to generate
Main:timesToShow = 10 ! show how far along run is this many times
Main:timesAllowErrors = 3 ! abort run after this many flawed events
Main:showChangedSettings = on ! print changed flags/modes/parameters
#Main:showAllSettings = on ! print all flags/modes/parameters
Main:showChangedParticleData = on ! print changed particle and decay data
#Main:showAllParticleData = on ! print all particle and decay data
 
! 2) Beam parameter settings. Values below agree with default ones.
Beams:idA = 2212 ! first beam, p = 2212, pbar = -2212
Beams:idB = 2212 ! second beam, p = 2212, pbar = -2212
Beams:eCM = 14000. ! CM energy of collision
 
! 3) Pick processes and kinematics cuts.
Charmonium:gg2QQbar[3S1(1)]g = on ! colour singlet charmonium production
Charmonium:gg2QQbar[3P0(1)]g = on ! "
Charmonium:gg2QQbar[3P1(1)]g = on ! "
Charmonium:gg2QQbar[3P2(1)]g = on ! "
PhaseSpace:pTHatMin = 20. ! minimum pT of hard process
 
! 4) Alternative beam and process selection from a Les Houches Event File.
! NOTE: to use this option, comment out the lines in section 3 above
! and uncomment the ones below. Section 2 is ignored for frameType = 4.
#Beams:frameType = 4 ! read info from a LHEF
#Beams:LHEF = ttbar.lhe ! the LHEF to read from
New file
/trunk/examples/pythia8/main31.cc
@@ -0,0 +1,92 @@
// main31.cc is a part of the PYTHIA event generator.
// Copyright (C) 2011 Mikhail Kirsanov, Torbjorn Sjostrand.
// PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
// Please respect the MCnet Guidelines, see GUIDELINES for details.
 
// This is a simple test program.
// It illustrates how HepMC can be interfaced to Pythia8.
// It studies the charged multiplicity distribution at the LHC.
// HepMC events are output to the hepmcout31.dat file.
// Written by Mikhail Kirsanov based on main01.cc.
 
// WARNING: typically one needs 25 MB/100 events at the LHC.
// Therefore large event samples may be impractical.
 
#include "Pythia.h"
#include "HepMCInterface.h"
 
#include "HepMC/GenEvent.h"
#include "HepMC/IO_GenEvent.h"
 
// Following line is a deprecated alternative, removed in recent versions
//#include "HepMC/IO_Ascii.h"
//#include "HepMC/IO_AsciiParticles.h"
 
// Following line to be used with HepMC 2.04 onwards.
#ifdef HEPMC_HAS_UNITS
#include "HepMC/Units.h"
#endif
 
using namespace Pythia8;
 
int main() {
 
// Interface for conversion from Pythia8::Event to HepMC one.
HepMC::I_Pythia8 ToHepMC;
// ToHepMC.set_crash_on_problem();
 
// Specify file where HepMC events will be stored.
HepMC::IO_GenEvent ascii_io("hepmcout31.dat", std::ios::out);
// Following line is a deprecated alternative, removed in recent versions
// HepMC::IO_Ascii ascii_io("hepmcout31.dat", std::ios::out);
// Line below is an eye-readable one-way output, uncomment the include above
// HepMC::IO_AsciiParticles ascii_io("hepmcout31.dat", std::ios::out);
 
// Generator. Process selection. LHC initialization. Histogram.
Pythia pythia;
pythia.readString("HardQCD:all = on");
pythia.readString("PhaseSpace:pTHatMin = 20.");
pythia.init( 2212, 2212, 14000.);
Hist mult("charged multiplicity", 100, -0.5, 799.5);
 
// Begin event loop. Generate event. Skip if error. List first one.
for (int iEvent = 0; iEvent < 100; ++iEvent) {
if (!pythia.next()) continue;
if (iEvent < 1) {pythia.info.list(); pythia.event.list();}
 
// Find number of all final charged particles and fill histogram.
int nCharged = 0;
for (int i = 0; i < pythia.event.size(); ++i)
if (pythia.event[i].isFinal() && pythia.event[i].isCharged())
++nCharged;
mult.fill( nCharged );
 
// Construct new empty HepMC event.
#ifdef HEPMC_HAS_UNITS
// This form with arguments is only meaningful for HepMC 2.04 onwards,
// and even then unnecessary if HepMC was built with GeV and mm as units..
HepMC::GenEvent* hepmcevt = new HepMC::GenEvent(
HepMC::Units::GEV, HepMC::Units::MM);
#else
// This form is needed for backwards compatibility.
// In HepMCInterface.cc a conversion from GeV to MeV will be done.
HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
#endif
 
// Fill HepMC event, including PDF info.
ToHepMC.fill_next_event( pythia, hepmcevt );
// This alternative older method fills event, without PDF info.
// ToHepMC.fill_next_event( pythia.event, hepmcevt );
 
// Write the HepMC event to file. Done with it.
ascii_io << hepmcevt;
delete hepmcevt;
 
// End of event loop. Statistics. Histogram.
}
pythia.statistics();
cout << mult;
 
// Done.
return 0;
}
New file
/trunk/examples/pythia8/CMakeLists.txt
@@ -0,0 +1,20 @@
 
set( example_code
main31.cc
main32.cc
main32.cmnd
README
)
 
CONFIGURE_FILE( ${CMAKE_CURRENT_SOURCE_DIR}/GNUmakefile.example.in
${CMAKE_CURRENT_BINARY_DIR}/GNUmakefile @ONLY )
CONFIGURE_FILE( ${CMAKE_CURRENT_SOURCE_DIR}/config.csh.in
${CMAKE_CURRENT_BINARY_DIR}/config.csh @ONLY )
CONFIGURE_FILE( ${CMAKE_CURRENT_SOURCE_DIR}/config.sh.in
${CMAKE_CURRENT_BINARY_DIR}/config.sh @ONLY )
 
INSTALL (FILES ${example_code}
${CMAKE_CURRENT_BINARY_DIR}/GNUmakefile
${CMAKE_CURRENT_BINARY_DIR}/config.csh
${CMAKE_CURRENT_BINARY_DIR}/config.sh
DESTINATION share/HepMC/examples/pythia8 )
New file
/trunk/examples/pythia8/README
@@ -0,0 +1,18 @@
These examples are from pythia8
 
To build and run the examples:
define PYTHIA8_HOME
check config.csh, config.sh, and GNUmakefile
edit variables if necessary
source the config file
make
./main31.exe
./main32.exe main32.cmnd hepmcout32.dat
 
bash:
export PYTHIA8_HOME=/some/path/to/pythia8xxx
source config.sh
 
csh:
setenv PYTHIA8_HOME /some/path/to/pythia8xxx
source config.csh
New file
/trunk/examples/pythia8/main32.cc
@@ -0,0 +1,141 @@
// main32.cc is a part of the PYTHIA event generator.
// Copyright (C) 2011 Mikhail Kirsanov, Torbjorn Sjostrand.
// PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
// Please respect the MCnet Guidelines, see GUIDELINES for details.
 
// This is a simple test program.
// It illustrates how a file with HepMC events can be generated by Pythia8.
// Input and output files are specified on the command line, e.g. like
// ./main32.exe main32.cmnd hepmcout32.dat > out
// The main program contains no analysis; this is intended to happen later.
// It therefore "never" has to be recompiled to handle different tasks.
 
// WARNING: typically one needs 25 MB/100 events at the LHC.
// Therefore large event samples may be impractical.
 
#include "Pythia.h"
#include "HepMCInterface.h"
 
#include "HepMC/GenEvent.h"
#include "HepMC/IO_GenEvent.h"
 
// Following line is a deprecated alternative, removed in recent versions.
//#include "HepMC/IO_Ascii.h"
//#include "HepMC/IO_AsciiParticles.h"
 
// Following line to be used with HepMC 2.04 onwards.
#ifdef HEPMC_HAS_UNITS
#include "HepMC/Units.h"
#endif
 
using namespace Pythia8;
 
int main(int argc, char* argv[]) {
 
// Check that correct number of command-line arguments
if (argc != 3) {
cerr << " Unexpected number of command-line arguments. \n You are"
<< " expected to provide one input and one output file name. \n"
<< " Program stopped! " << endl;
return 1;
}
 
// Check that the provided input name corresponds to an existing file.
ifstream is(argv[1]);
if (!is) {
cerr << " Command-line file " << argv[1] << " was not found. \n"
<< " Program stopped! " << endl;
return 1;
}
 
// Confirm that external files will be used for input and output.
cout << " PYTHIA settings will be read from file " << argv[1] << endl;
cout << " HepMC events will be written to file " << argv[2] << endl;
 
// Interface for conversion from Pythia8::Event to HepMC one.
HepMC::I_Pythia8 ToHepMC;
// ToHepMC.set_crash_on_problem();
 
// Specify file where HepMC events will be stored.
HepMC::IO_GenEvent ascii_io(argv[2], std::ios::out);
// Following line is a deprecated alternative, removed in recent versions
// HepMC::IO_Ascii ascii_io("hepmcout32.dat", std::ios::out);
// Line below is an eye-readable one-way output, uncomment the include above
// HepMC::IO_AsciiParticles ascii_io("hepmcout32.dat", std::ios::out);
// Generator.
Pythia pythia;
 
// Read in commands from external file.
pythia.readFile(argv[1]);
 
// Extract settings to be used in the main program.
int nEvent = pythia.mode("Main:numberOfEvents");
int nShow = pythia.mode("Main:timesToShow");
int nAbort = pythia.mode("Main:timesAllowErrors");
bool showCS = pythia.flag("Main:showChangedSettings");
bool showAS = pythia.flag("Main:showAllSettings");
bool showCPD = pythia.flag("Main:showChangedParticleData");
bool showAPD = pythia.flag("Main:showAllParticleData");
// Initialization. Beam parameters set in .cmnd file.
pythia.init();
 
// List settings.
if (showCS) pythia.settings.listChanged();
if (showAS) pythia.settings.listAll();
 
// List particle data.
if (showCPD) pythia.particleData.listChanged();
if (showAPD) pythia.particleData.listAll();
 
// Begin event loop.
int nPace = max(1, nEvent / max(1, nShow) );
int iAbort = 0;
for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
if (nShow > 0 && iEvent%nPace == 0)
cout << " Now begin event " << iEvent << endl;
 
// Generate event.
if (!pythia.next()) {
 
// If failure because reached end of file then exit event loop.
if (pythia.info.atEndOfFile()) {
cout << " Aborted since reached end of Les Houches Event File\n";
break;
}
 
// First few failures write off as "acceptable" errors, then quit.
if (++iAbort < nAbort) continue;
cout << " Event generation aborted prematurely, owing to error!\n";
break;
}
 
// Construct new empty HepMC event.
#ifdef HEPMC_HAS_UNITS
// This form with arguments is only meaningful for HepMC 2.04 onwards,
// and even then unnecessary if HepMC was built with GeV and mm as units..
HepMC::GenEvent* hepmcevt = new HepMC::GenEvent(
HepMC::Units::GEV, HepMC::Units::MM);
#else
// This form is needed for backwards compatibility.
// In HepMCInterface.cc a conversion from GeV to MeV will be done.
HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
#endif
 
// Fill HepMC event, including PDF info.
ToHepMC.fill_next_event( pythia, hepmcevt );
// This alternative older method fills event, without PDF info.
// ToHepMC.fill_next_event( pythia.event, hepmcevt );
 
// Write the HepMC event to file. Done with it.
ascii_io << hepmcevt;
delete hepmcevt;
 
// End of event loop. Statistics.
}
pythia.statistics();
 
// Done.
return 0;
}
/trunk/examples/CMakeLists.txt
@@ -16,3 +16,4 @@
DESTINATION share/HepMC/examples )
 
add_subdirectory(fio)
add_subdirectory(pythia8)