File: pythia8_example.cc

package info (click to toggle)
hepmc3 3.1.2-2.1
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 14,124 kB
  • sloc: fortran: 66,849; cpp: 13,604; ansic: 1,374; xml: 109; sh: 72; makefile: 33
file content (61 lines) | stat: -rw-r--r-- 1,807 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
/**
 *  @example pythia8_example.cc
 *  @brief Basic example of use for pythia8 interface
 *
 */
#include "HepMC3/GenEvent.h"
#include "HepMC3/WriterAscii.h"
#include "HepMC3/Print.h"

#include "Pythia8/Pythia.h"
#include "Pythia8ToHepMC3.h"

#include <iostream>
using namespace HepMC3;

/** Main program */
int main(int argc, char **argv) {
    if (argc < 3) {
        cout << "Usage: " << argv[0] << " <pythia_config_file> <output_hepmc3_file>" << endl;
        exit(-1);
    }

    Pythia8::Pythia pythia;
    Pythia8ToHepMC3 pythiaToHepMC;
    pythia.readFile(argv[1]);
    pythia.init();
    shared_ptr<GenRunInfo> run = make_shared<GenRunInfo>();
    struct GenRunInfo::ToolInfo generator={std::string("Pythia8"),std::to_string(PYTHIA_VERSION).substr(0,5),std::string("Used generator")};
    run->tools().push_back(generator);
    struct GenRunInfo::ToolInfo config={std::string(argv[1]),"1.0",std::string("Control cards")};
    run->tools().push_back(config);
    std::vector<std::string> names;
    for (int iWeight=0; iWeight < pythia.info.nWeights(); ++iWeight) {
     std::string s=pythia.info.weightLabel(iWeight);
     if (!s.length()) s=std::to_string((long long int)iWeight);
     names.push_back(s);
    }
    if (!names.size()) names.push_back("default");
    run->set_weight_names(names);
    WriterAscii file(argv[2],run);

    int nEvent = pythia.mode("Main:numberOfEvents");

    for( int i = 0; i< nEvent; ++i ) {
        if( !pythia.next() ) continue;

        GenEvent hepmc( Units::GEV, Units::MM );

        pythiaToHepMC.fill_next_event(pythia.event, &hepmc, -1, &pythia.info);

        if( i==0 ) {
            std::cout << "First event: " << std::endl;
            Print::listing(hepmc);
        }

        file.write_event(hepmc);
    }

    file.close();
    pythia.stat();
}