This file is indexed.

/usr/share/doc/pythia8-doc/html/HepMCInterface.html is in pythia8-doc-html 8.1.80-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
<html>
<head>
<title>HepMC Interface</title>
<link rel="stylesheet" type="text/css" href="pythia.css"/>
<link rel="shortcut icon" href="pythia32.gif"/>
</head>
<body>

<h2>HepMC Interface</h2>

An interface to the HepMC [<a href="Bibliography.html" target="page">Dob01</a>] standard event record 
format has been provided by M. Kirsanov. To use it, the relevant 
libraries need to be linked, as explained in the <code>README</code> 
and <code>README.HepMC</code> files. Only version 2.06 (or later) 
of HepMC is supported as of 1 January 2013, by agreement with the 
LHC experimental community.

<p/>
The (simple) procedure to translate PYTHIA 8 events into HepMC ones 
is illustrated in the <code>main41.cc</code>, <code>main42.cc</code>
<code>main61.cc</code> and <code>main62.cc</code>   
main programs. At the core is a call to the
<pre>
HepMC::Pythia8ToHepMC::fill_next_event( pythia, hepmcevt, ievnum = -1) 
</pre>
which takes a reference of the generator object and uses it, on the one
hand, to read out and convert the event record in <code>pythia.event</code> 
and, on the other hand, to extract and store parton-density (PDF),
cross section and other information for the hard subprocess from 
<code>pythia.info</code>. There is also an alternative form that
does not requires access to the full <code>pythia</code> object, 
but only the event record, at the expense of a reduced information
storage, see below.

<p/>
While PYTHIA always stores momenta in units of GeV, with <i>c = 1</i>, 
HepMC nowadays can be built either for MeV or GeV as default, a choice
that can then be overridden on an event-by-event basis, see e.g. the 
<code>main41.cc</code> code. When filling the HepMC event record, PYTHIA 
will convert to the unit specified for the current HepMC event record. 

<p/>
Also for length units there are choices, and again the PYTHIA interface 
will convert to the units set for the HepMC event record. Here the mm
choice of PYTHIA seems to be shared by most other programs, and is
HepMC default. 

<p/>
The status code is now based on the new standard introduced for HepMC 2.05,
see the <a href="EventRecord.html" target="page">Event::statusHepMC(...)</a> 
conversion routine for details. 

<p/>
The current values from <code>pythia.info.sigmaGen()</code> and 
<code>pythia.info.sigmaErr()</code> are stored for each event,
multiplied by <i>10^9</i> to convert from mb to pb. Note that
PYTHIA improves its accuracy by Monte Carlo integration in the course
of the run, so the values associated with the last generated event
should be the most accurate ones. If events also come with a dimensional 
weight, like in some Les Houches strategies, this weight is in units of pb.

<h2>The public methods</h2>

Here comes a complete list of all public methods of the 
<code>Pythia8ToHepMC</code> class in the <code>HepMC</code> 
(<i>not</i> <code>Pythia8</code>!) namespace.

<a name="method1"></a>
<p/><strong>Pythia8ToHepMC::Pythia8ToHepMC() &nbsp;</strong> <br/>
  
<strong>virtual Pythia8ToHepMC::~Pythia8ToHepMC() &nbsp;</strong> <br/>
the constructor and destructor take no arguments.
  

<a name="method2"></a>
<p/><strong>bool Pythia8ToHepMC::fill_next_event( Pythia8::Pythia& pythia, GenEvent* evt, int ievnum = -1) &nbsp;</strong> <br/>
convert a <code>Pythia</code> event to a <code>HepMC</code> one.
Will return true if succeeded.
<br/><code>argument</code><strong> pythia </strong>  : 
the input <code>Pythia</code> generator object, from which both the 
event and other information can be obtained.
   
<br/><code>argument</code><strong> evt </strong>  : 
the output <code>GenEvt</code> event, in its standard form.
   
<br/><code>argument</code><strong> iev </strong>  : 
set the event number of the current event. If negative then the 
internal event number is used, which is incremented by one for
each new event.
  
  

<a name="method3"></a>
<p/><strong>bool Pythia8ToHepMC::fill_next_event( Pythia8::Event& pyev, GenEvent* evt, int ievnum = -1, Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0) &nbsp;</strong> <br/>
convert a <code>Pythia</code> event to a <code>HepMC</code> one.
Will return true if succeeded.
<br/><code>argument</code><strong> pyev </strong>  : 
the input <code>Pythia</code> event that is to be converted to HepMC
format.
   
<br/><code>argument</code><strong> evt </strong>  : 
the output <code>GenEvt</code> event, in its standard form.
   
<br/><code>argument</code><strong> iev </strong>  : 
set the event number of the current event. If negative then the 
internal event number is used, which is incremented by one for
each new event.
  
<br/><code>argument</code><strong> pyinfo </strong>  : 
pointer to the <code>Pythia Info</code> object, which is used to 
extract PFD values, and process and cross section information.
Without such a pointer this information therefore cannot be stored,
i.e. it is equivalent to the three <code>set_store</code> methods
below being set false. 
   
<br/><code>argument</code><strong> pyset </strong>  : 
pointer to the <code>Pythia Settings</code> object, which is used to 
decide whether hadronization is carried out, and therefore whether
to warn about unhadronized partons. Without such a pointer the
<code>set_free_parton_warnings</code> method below uniquely controls
the behaviour.
   
  

<p/>
The following paired methods can be used to set and get the status of 
some switches that modify the behaviour of the conversion routine. 
The <code>set</code> methods have the same default input values as 
the internal initialization ones, so that a call without an argument 
(re)stores the default.

<a name="method4"></a>
<p/><strong>void Pythia8ToHepMC::set_print_inconsistency(bool b = true) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::print_inconsistency() &nbsp;</strong> <br/>
print a warning line, on <code>cerr</code>, when inconsistent mother 
and daughter information is encountered.
  

<a name="method5"></a>
<p/><strong>void Pythia8ToHepMC::set_free_parton_warnings(bool b = true) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::free_parton_warnings() &nbsp;</strong> <br/>
check and print a warning line when unhadronized gluons or quarks are 
encountered in the event record. Does not apply when Pythia hadronization 
is switched off. Default is to do this check.
  

<a name="method6"></a>
<p/><strong>void Pythia8ToHepMC::set_crash_on_problem(bool b = false) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::crash_on_problem() &nbsp;</strong> <br/>
if problems (like the free partons above) are encountered then the run 
is interrupted by an <code>exit(1)</code> command. Default is not to crash.
  

<a name="method7"></a>
<p/><strong>void Pythia8ToHepMC::set_convert_gluon_to_0(bool b = false) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::convert_gluon_to_0() &nbsp;</strong> <br/>
the normal gluon identity code 21 is used also when parton density
information is stored, unless this optional argument is set true to
have gluons represented by a 0. This choice does not affect the 
normal event record, where a gluon is always 21. 
  

<a name="method8"></a>
<p/><strong>void Pythia8ToHepMC::set_store_pdf(bool b = true) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::store_pdf() &nbsp;</strong> <br/>
for each event store information on the two incoming flavours, their
x values and common factorization scale, and the values of the two 
parton distributions, <i>xf(x,Q)</i>.
  

<a name="method9"></a>
<p/><strong>void Pythia8ToHepMC::set_store_proc(bool b = true) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::store_proc() &nbsp;</strong> <br/>
for each event store information on the Pythia process code, the 
renormalization scale, and <i>alpha_em</i> and <i>alpha_s</i>
values used for the hard process. 
  

<a name="method10"></a>
<p/><strong>void Pythia8ToHepMC::set_store_xsec(bool b = true) &nbsp;</strong> <br/>
  
<strong>bool Pythia8ToHepMC::store_xsec() &nbsp;</strong> <br/>
for each event store information on the Pythia cross section and its error,
in pb, and the event weight. If events also come with a dimensional weight, 
like in some Les Houches strategies, this weight is in units of pb.
  

</body>
</html>

<!-- Copyright (C) 2013 Torbjorn Sjostrand -->