| 5 |
Advanced usage of WHIZARD |
|
In the preceding chapters we have described the use of WHIZARD as a
stand-alone program. As such, it is able to calculate total cross
sections and to generate events on the partonic level, or even on the
hadronic level if fragmentation is enabled. The events can be
analyzed online (on the partonic level) or written to file (on the
partonic or on the hadronic level).
This already covers many real-life applications of a Monte-Carlo
integration and event-generation package. For the theorist, it may be
sufficient to calculate a total cross section with a given set of
parameters and to plot differential distributions. For an
experimental simulation, the generated event files can be used for
further processing.
However, one encounters situations where a more elaborated framework
is desired. For instance, one would like to consider multiple sets of
cuts, energies and parameter sets. For each set, the input file has
to be modified, the program has to be run, and the results have to be
stored somewhere in order to be compared later. This can easily
become cumbersome, and it is desirable to have a generic approach.
Another common problem is the selection of a certain subclass of
Feynman diagrams contributing to a given process, either for speeding
up the calculation by dropping irrelevant terms, or for getting a
better understanding of the relative weight of interfering
subprocesses. Some possibilities are described in the following
sections.
| 5.1 |
Feynman diagram selections |
|
While the matrix element generators called by WHIZARD, by default,
generate an amplitude which corresponds to the complete set of Feynman
diagrams which connect a given initial and final state, there are
problems where only a subset is needed. Therefore, WHIZARD has
support for diagram selections. Since O'Mega does not use Feynman
diagrams, this option is currently available only with the CompHEP and
MadGraph matrix element generators. (O'Mega supports
restrictions on the internal lines of an amplitude. This is explained
below in Sec. 5.2.)
The selection must be made before the final executable is generated.
In the configuration file, the process must be marked accordingly by
the letter ``s'' in the last column, for instance:
# Processes
# Tag In Out Method Option
#=====================================================
cc3 e1,E1 e2,N2,u,D chep s
Note that this makes sense only for the methods chep and
mad, as mentioned before. Next, the diagram numbers to be
selected have to be determined. To this end, for each process with
this option a diagram selection file has to be put into the
configuration directory.
Since the user usually does not know diagram numbers in advance, it is
possible to do an empty run of the matrix element generators first.
If the command
make diagrams
is issued, matrix element generation is started as usual, but the
Fortran95 source code is not actually generated (at least, with
CompHEP: for the other generators it makes little difference).
Instead, after this command is completed, you will find PostScript
files with Feynman diagrams in graphical representation in the
subdirectory process-src. View or print them and note, for
each process, the diagrams you are interested in. In the subdirectory
conf you should now find, for each process, a file with the
extension .sel, which contains one entry for each diagram.
Initially, all entries are set 0 (false). As a consequence,
if the file is not modified, a second WHIZARD run would fail since
all diagrams are deleted initially. Change all entries
you wish to enable into T or 1 (true). If you now
rerun WHIZARD to generate the actual source code
make proc
or even
make prg install
WHIZARD should read the selection files and take into account only
the specified subsets of diagrams. You may control this by looking at
the PostScript output this second pass will produce: it should contain
only the Feynman diagrams you have selected.
Needless to say, this feature should be used with great care. By
violating electroweak or even electromagnetic gauge invariance, the
results may be off by far, since delicate cancellations may have been
spoiled. As a rule, never delete diagrams in any gauge other than
unitarity gauge (in the CompHEP) case, and remember which diagrams
have been selected (i.e., save the PostScript file somewhere) when
storing results, since the output of the WHIZARD executable itself
does not contain any hint that only a subset of Feynman diagrams have
been taken into account. In many cases the same effect can also be
achieved in a simpler way by putting some parameter to extreme values.
For instance, instead of removing Higgs diagrams by hand, one may set
the Higgs mass to a very large value --- the results will then be
valid in any gauge.
| 5.2 |
Restrictions on intermediate states |
|
As an alternative (and often more convenient) method for selecting
Feynman diagrams, the O'Mega interface supports restrictions on
internal lines.18 The WHIZARD implementation of this interface is
experimental and may be improved or changed later. As an example,
the following specification selects the resonant WW diagrams,
analogous to the previous example:
# Processes
# Tag In Out Method Option
#=====================================================
cc3 e1,E1 e1,N1,u,D omega r: 3+4~W- && 5+6~W+
Simultaneous requirements for a process are concatenated by &&
(AND).
The following examples select only diagrams with an s-channel and
t-channel g*, respectively:
cc_gamma_s e1,E1 e1,N1,u,D omega r: 3+4+5+6~A
cc_gamma_t e1,E1 e1,N1,u,D omega r: 1+3~A
The in- and out-particles are numbered consecutively. The tilde
selects only diagrams where the specified combination of momenta
corresponds to the propagator for the given particle. A question mark
instead of a particle name stands for all possible propagators with
the specified momentum combination. These particle names have
to follow the O'Mega notation; they are not translated by WHIZARD
as in the process specification. In particular, for colored
amplitudes ('c' option) color indices have to be appended to particle
names, which makes the interface inconvenient for this case. It is
important to keep in mind that such a selection typically violates
gauge invariance, since the intermediate particle is not forced on
shell. (An on-shell option is also present in O'Mega, but not yet
supported by WHIZARD.)
The || (OR) relation is also possible, and sometimes
needed, as in this example:
ww_or_zz e1,E1 e1,N1,E1,n1 omega r: 3+4~W- && 5+6~W+ || 3+5~Z && 4+6~Z
AND and OR relations can be grouped by parantheses ().
Currently, WHIZARD does not interpret or use such restrictions for
its phase space setup, so the possible speedup is constrained to the
matrix element evaluation. Nevertheless, this feature can be very useful.
| 5.3 |
User-defined spectra and structure functions |
|
It is possible to insert arbitrary beam spectra into the WHIZARD
code. If the corresponding code is included in the user module,
it can be called from the main program. To do this, in the
whizard-src subdirectory, copy the file user.f90.default
to user.f90 and edit the appropriate functions defined within
this module according to your needs. The code you have written will
be compiled instead of the default code.
The comments in the user.f90.default source should guide you.
There is the possibility to define either single-beam (factorized) or
correlated spectra. The particle codes which have been specified in
the beam_input block of the input file can be accessed in order to
select between different spectra for different beams.
In the input file, for each beam there is a master switch
USER_spectrum_on which determines whether the user code is
called. If no user code has been written, by default a uniform
x spectrum between 0 and 1 is inserted which could be used for
testing purposes. Furthermore, two parameters can be specified:
USER_spectrum_mode, an integer which the user function can access
to select between different spectra for the same particle, and
USER_spectrum_sqrts which by default is the total c.m. energy,
but in principle can be used to transfer any real constant parameter.
Apart from user-defined spectra, user-defined structure functions can
be inserted. They are handled in an exactly analogous way (with
spectrum replaced by strfun in parameter names). The
only difference is that user spectra are applied before any
built-in spectra and structure functions (like beamstrahlung, PDFs,
etc.), while user structure functions are applied after built
in ones, i.e., they are called for the partonic initial state. For
applications, think of a photon collider spectrum for the initial
beams as opposed to the scattering of W bosons emitted from quarks.
The efficiency of the program in the presence of user-defined spectra
or structure functions depends on the way they are implemented by the
user. It is straightforward to insert the corresponding functions in
a literal way as functions f(x) where x is the beam energy
fraction. However, if f is strongly peaked this might be lead to
poor performance. In many cases one rather should apply a variable
transformation
|
|
ó
õ |
|
dx f(x) s(xs) =
|
ó
õ |
|
dy |
|
f(g(y)) s(xs)
(8) |
where x=g(y). The user interface allows for this, where the input
parameters are identified with y and the output parameters x=g(y)
can be different. In that case, the spectral function f(x)=f(g(y))
has to be multiplied by the Jacobian dg/dy within the structure
function code.
The implementation can be checked by making use of the ``test'' matrix
elements (Sec. 4.1.2). You may set up a placeholder process
like
sftest e1,E1 A,A test
where the final state is massless and remove all cuts (in the example,
set default_Q_cut=0 or define an empty entry in the cut file).
The matrix element is set constant, and the would-be cross section for
the process is the integral of the structure function
which is usually known. Furthermore, the energy distribution of the
final state should reproduce the structure function shape.
In the current implementation, beam polarization (if any) cannot be
accessed or modified by the user spectra or structure functions. Beam
polarization will be kept unchanged if the in- and out-particle of the
corresponding splitting coincide (as for ISR), and it will be ignored
in the opposite case (as for PDFs).
In many applications the built-in capabilities for cuts are not
sufficient. If cuts are defined using the standard interface, all
cuts are composed such that an event is only retained if it passes all
cuts (logical and). There can be several windows in one cut
variable, but beyond that it is not possible to compose cuts by a
logical or. Furthermore, you may be missing your favorite cut
model in the list of predefined variables.
User-defined cut algorithms can easily be placed in the user
module which comes with the distribution, the same module where
user-defined spectra can be coded. To activate the extra code, the
parameter user_cut_mode in the input file should be set to a
nonzero integer. Then, the new cut function will be applied
throughout program execution in addition to the cut conditions defined
in the standard way. Within the cut function, the numerical value of
user_cut_mode is used to select one of several cut modes
defined by the user. (To disable the default cuts, which may be
desired in this case, write an empty process entry in
whizard.cut1.)
The cut function takes an array of momenta (this type and the
functions operating on it is defined in the module lorentz.f90)
and results in the logical variable accept. The default behavior
is to unconditionally set accept=.true., but user code would
change that, deciding whether to accept the event based on the
momentum configuration and the particle codes.
As an example, let us implement a cut such that an event is accepted
if there is a photon with energy larger than 50 GeV. (This is not
possible by the built-in cut capabilities, since for a generic process
we do not know which photon to check.) To this end, we first make a
copy of the standard (no-op) user file, which then can be changed:
cd whizard-src
cp user.f90.default user.f90
The file user.f90 then has to be edited. We take the section
which refers to the user cut definition, which originally reads
! This minimal cut function accepts everything
function cut (p, code, mode) result (accept)
type(four_momentum), dimension(:), intent(in) :: p
integer, dimension(:), intent(in) :: code
integer, intent(in) :: mode
logical :: accept
select case (mode)
case default
accept = .true.
end select
end function cut
and change it as follows:
! This cut function accepts an event only if there is a photon with
! energy larger than 50 GeV
function cut (p, code, mode) result (accept)
type(four_momentum), dimension(:), intent(in) :: p
integer, dimension(:), intent(in) :: code
integer, intent(in) :: mode
logical :: accept
integer :: i
select case (mode)
case (1)
accept = .false.
do i = 1, size (code)
if (code(i) == PHOTON .and. energy (p(i)) > 50) then
accept = .true.
end if
end do
case default
accept = .true.
end select
end function cut
We have inserted a loop over all outgoing particles, such that the
code can be applied to an arbitrary process. For each particle
i, its PDG code is checked against the predefined code
PHOTON (which is defined as 22 in particles.f90), and
the function energy is evaluated for the corresponding momentum
p(i). This function, which takes a variable of type
four_momentum as argument and returns a real number of default
kind (i.e., double precision), is defined in lorentz.f90.
When WHIZARD is now (re-)compiled and installed,
make install
the user code becomes part of the program. To activate it, the
integration_input in the whizard.in file has to specify
a nonzero user-cut mode, i.e.,
&integration_input
...
user_cut_mode = 1
/
since in the function defined above, the mode 1 refers to
the special treatment of photons. When the program is started with
this setting, you will see the message
! Activating additional preselection cuts in user.f90: Mode # 1
in the screen output, and the cross section and event distribution
should reflect the new cut condition.
| 5.5 |
Reweighting matrix elements |
|
Beyond the possibility to apply arbitrary cuts, the user might be
interested in reweighting events. One should distinguish two different
situations:
-
A given event sample is rescanned, recalculating the matrix element
values event by event with some modification applied. The resulting
events carry a weight which is given by the ratio of the new matrix
element values compared to the old ones. The event selection is still
based on the unmodified values.
To achieve this, one should first generate a reference sample in the
usual way (adaptation and event generation) which results in the raw
event file whizard.evx. Then, one may change the input
parameters as desired, set the parameter recalculate=T, and
start a second run. This will have two effects: First, the parameters
read_grids and read_events are assumed also to be set
(even if the user has not modified them), together with
read_grids_force and read_events_force.19 Thus, the grids and events of the
first run are reused. Second, when the event sample is read, the
particle momenta are left unchanged, but the matrix element values and
the quantum number probabilities are recalculated, using the new input
data.
The change in the matrix element values will be taken into account in
the internal analysis: The integral and histograms will be calculated
using the ratio of new and old matrix element values as event
weights. If the results of the reference run have been stored (one
should rename output files to prevent them from being overwritten),
the results can easily be compared.
To make further use of the reweighted events, one can write the events
to file (write_events=T) in ``long'' (ASCII) format:
write_events_format=3. In addition to the particle codes,
helicities and momenta, the file will contain the (recalculated)
matrix element value and the ratio of new and old values for each
event. Again, if several parameter sets are to be compared, one
should make sure to rename the output file whizard.evt after
each run to prevent overwriting.
In any case, one has to be careful. While it does not harm to change
a coupling constant by a small amount (e.g., an anomalous coupling)
and rescan the events, modifying particle masses and cuts may have
undesired side-effects since the momentum configurations are left
intact. For instance, with new mass assignments particles may no
longer be on-shell, such that delicate gauge cancellations in the
matrix elements are invalid, resulting in unphysical effects.
- The matrix element function itself is modified by a user-defined
weight, a function of the outgoing momenta. The weight function is
applied throughout adaptation and event generation, so the integral
and the event selection depend on it.
A natural application is the calculation of integrated observables
other than the total cross section, e.g. moments of the angular
distribution. One might also think of radiative corrections in
``K-factor'' form: a function of the momenta which multiplies the
lowest-order matrix element.20 The corresponding
code can be inserted in the weight routine in the
user.f90 module. Similar to the cut routine
(Sec.5.4), it takes the array of particle momenta and
particle codes to generate a single number weight. In the
default routine, weight is always equal to one. The user may
change this to an arbitrary function, provided it is nonnegative in
the whole accessible phase space.21 The
user-defined weight function is then activated by setting
user_weight_mode nonzero in the input file. The numerical
value can be used to select between different reweighting modes.
It is perfectly possible to combine the two approaches: One may define
one or more reweighting functions in the user module. As a
reference run, events are generated with user_weight_mode=0.
Then, both recalculate is set and the mode is set nonzero in a
second run. A user-defined cut function can be treated
analogously.
WHIZARD allows to specify options on the command line.22 Options are available both in long form (double
dashes) and in short form (a single dash). For instance, the
following command executes WHIZARD in subdirectory tmp for
the process ee_nnh, without reading an input file:
./whizard -nd tmp --process_input 'process_id = "ee_nnh"'
These options are executed by WHIZARD before any files are
read or written
|
| Option |
Long version |
Value |
Description |
|
-h |
--help |
|
Display usage message and exit. |
-V |
--version |
|
Display version string and exit. |
-l |
--list |
|
List all available processes and exit. |
-d |
--directory |
directory |
Set working
directory. |
-f |
--filename |
file name |
Set default base name
for all input and output files (without extension). |
-i |
--input_file |
file name |
Set base name for
the main input file (without extension). This overrides the -f option. |
-n |
--no_input_file |
|
Do not read the main input file. |
|
If a working directory different from the default one is chosen by the
-d option, one has to make sure that either a phase space file
(.phs or .phx) or the model definition file
(.mdl) is present in the chosen directory. Otherwise, one
could use command-line options to point to files in the initial
directory (or anywhere else), for instance:
./whizard --directory tmp --integration_input 'read_model_file="./whizard"'
When a --directory setting is in effect, it can be overriden by
filenames that begin with ``/'' (absolute file name) or with
``./'' (filename relative to the initial directory).
Therefore, in the example the model file whizard.mdl will be
looked for in the directory where the command was executed, while all
other files will be looked for in or written to the subdirectory
tmp.
| 5.6.2 |
Options for setting parameters |
|
Parameters can be set not just in the input file, but also on the
command line. The command-line parameters are read after all
input files have been read; they override settings in the
input file(s). Thus, one can set up a generic input file
whizard.in where the constant parameters are defined (or no
input file at all) and use the command line to run WHIZARD for,
e.g., a specific process, beam energy and Higgs boson mass:
./whizard -p 'process_id="ee_nnh" sqrts=800' -P 'mH = 130'
Note that the option arguments are interpreted verbatim as if they
were contained in the input file, so string values have to be enclosed
by extra quotes.
|
| Option |
Long version |
Value |
Description |
|
-p |
--process_input |
parameters |
set
parameters in the &process_input block. |
-I |
--integration_input |
parameters |
set
parameters in the &integration_input block. |
-S |
--simulation_input |
parameters |
set
parameters in the &simulation_input block. |
-D |
--diagnostics_input |
parameters |
set
parameters in the &diagnostics_input block. |
-P |
--parameter_input |
parameters |
set
parameters in the ¶meter_input block. |
-B1 |
--beam_input1 |
parameters |
set
parameters in the first &beam_input block. |
-B2 |
--beam_input2 |
parameters |
set
parameters in the second &beam_input block. |
|
| 5.7 |
Tables and file management |
|
The WHIZARD program, after it has been compiled with a specific list
of processes, may be viewed as a kind of filter: It takes
as input
-
the model definition file
whizard.mdl,
- the parameter definition file
whizard.in, and optionally
- the cut definition file
whizard.cut1
- and the cut/analysis definition file
whizard.cut5,
and transforms the information contained within into output files,
namely
-
the integration results and diagnostics, collected in the file
whizard.out,
- the integration grids, usable to bypass further time-consuming
adaptation runs (
whizard.grb), and optionally
- an event file usable for re-running the program (
whizard.evx);
- an event file for further processing by external programs
(
whizard.evt);
- and a set of histograms (
whizard-plots.dat) which can be
converted to PostScript format (whizard-plots.ps).
For a second WHIZARD run, the input files have to be edited, and the
output files will be overwritten if they have not been moved to a
different location. Instead of doing this by hand, one could write a
script that executes a loop where in each iteration the relevant line
in the input file is modified, the program is run, and the output
files are renamed or moved such that they do not get overwritten.
Here is such a script:
#!/usr/bin/perl
# Loop over the total energy, write a table and store the results
# in subdirectories
$sqrts_min = 350; $sqrts_max = 800; $sqrts_step = 50;
# These files will be preserved
@outfiles = ("whizard.out", "whizard.grb", "whizard.evt");
# This will become the table of results
$resultfile = "whizard.out";
$table = sprintf(" %-10s %-20s %-20s\n", "sqrts [GeV]",
"Total cross section [fb]", "Integration error [fb]");
# Main loop
for ($sqrts=$sqrts_min; $sqrts<=$sqrts_max; $sqrts=$sqrts+$sqrts_step) {
# Read the input file template, modify the sqrts entry, and write the
# new input file
$infile = "whizard.in.0"; $outfile = "whizard.in";
open(INFILE, $infile) or die "Can't open $infile for reading";
open(OUTFILE, ">$outfile") or die "Can't open $outfile for writing";
while (<INFILE>) {
s/sqrts = \d+/sqrts = $sqrts/;
print OUTFILE;
}
close(OUTFILE);
close(INFILE);
# Run WHIZARD
print "\n**** sqrts = $sqrts ****\n";
$stat = system("./whizard");
$stat==0 or die "WHIZARD run exited abnormally";
# Extract the integration results and add them to the table
open(RESULTS, "grep 'integral =' $resultfile|")
or die "Grep failed on $resultfile";
while(<RESULTS>) { ($dummy, $equals, $integral) = split; }
close(RESULTS);
open(RESULTS, "grep 'error =' $resultfile|")
or die "Grep failed on $resultfile";
while(<RESULTS>) { ($dummy, $equals, $error) = split; }
close(RESULTS);
$table .= sprintf(" %10g %20g %20g\n",
$sqrts, $integral, $error);
# Create a subdirectory and move the results there
$dirname = "sqrts=$sqrts";
mkdir($dirname, oct("0755")) or die "Can't mkdir $dirname";
foreach $file(@outfiles) {
system("mv", $file, $dirname);
}
}
print "\n*** Results:\n";
print $table;
print "*** End of WHIZARD runs\n";
exit 0;
The script assumes that the subdirectories do not yet exist, and that
you have a template whizard.in.0 for the input file which
contains an entry like sqrts = 0. The table is shown on screen
after all runs have been completed.
This kind of task can be formalized further. The program
whizwrap (written by Th. Ohl) takes care of loops over
parameters and provides a convenient management of output files. It
is called from the command line; the kind of task to be performed by
WHIZARD is specified in terms of command-line options. A future
version may provide a graphical interface.
| 5.8 |
WHIZARD as a subroutine library |
|
A more ambitious enterprise is the integration of WHIZARD into a
larger framework which may also account for the subsequent steps of
experimental simulation and analysis. Even if this is not desired, it
may be useful to access WHIZARD by subroutine calls in wrapper
programs written in Fortran, C, or other program languages.
To make this possible, WHIZARD is organized as a set of libraries
(in the UNIX sense). When the configuration file (the process list)
has been specified, instead of making a stand-alone executable by
make prg
one may generate the libraries only
make lib
which then can be linked with a user-supplied main program.
To facilitate this procedure, after executing make lib, the
user will find the generated libraries in the subdirectory lib
(except for system-supplied libraries from the CERNLIB), together with
a script called whizard.ld in the bin subdirectory which
links an arbitrary program with the libraries in correct order. Thus,
the user may write his own Fortran 95 main program, call it, e.g.,
my_main.f90, compile it23 and execute
bin/whizard.ld -o my_whizard my_main.o
to generate his own version of WHIZARD.24
Here is a main program that performs the same task as the script
shown in the preceding section:
! Loop over the total energy, write a table and store the results in
! separate output files.
program my_main
use kinds ! Defines the default real kind
use whizard ! Defines the WHIZARD subroutines
implicit none ! No implicit typing
! Parameters
integer, parameter :: n_values = 10
real(kind=default), parameter :: sqrts_min = 350, sqrts_max = 800
! Variables
integer :: i
character(len=14+3) :: filename
! Tables
real(kind=default), dimension(n_values) :: &
sqrts, integral, error, chi2, efficiency
! Main loop
do i = 1, n_values
! Calculate energy, set filename accordingly
sqrts(i) = ((i-1) * sqrts_max + (n_values-i) * sqrts_min) / (n_values-1)
write(filename, "(A14,I3)") "whizard_sqrts=", nint(sqrts(i))
call whizard_set_filename (filename)
! Read input and override energy setting
write (*,*)
write (*, "(A,1x,I3,1x,A)") "*** sqrts =", nint(sqrts(i)), "***"
call whizard_read_input
whizard_input%sqrts = sqrts(i)
! Integrate, record result, generate events
call whizard_integrate
call whizard_get_results (integral(i), error(i), chi2(i), efficiency(i))
call whizard_generate
call whizard_end
end do
! Write table
write (*,*)
write (*, "(A)") "*** Results:"
write (*, "(1x,A15,2x,A25,2x,A25)") "sqrts [GeV]", &
"Total cross section [fb]", "Integration error [fb]"
do i = 1, n_values
write (*, "(1x,G15.5,2x,G25.5,2x,G25.5)") &
sqrts(i), integral(i), error(i)
end do
write (*, "(A)") "*** End of WHIZARD runs"
end program my_main
The meaning of the individual subroutine calls is described in detail
in the following section. Note that parameters in the input block can
be accessed via the object whizard_input which is a module
variable of the whizard module. More exactly, variables in the
process_input, integration_input,
simulation_input and diagnostics_input blocks can be
accessed as above
whizard_input%sqrts = sqrts
while variables in the parameter_input block are accessed like
whizard_input%par%mh = 115
and beam parameters are accessed like
whizard_input%beam(1)%PDF_nset = 43
| 5.9 |
Interface of the WHIZARD module |
|
All subroutines described in this section are are accessed via the
module whizard, which is included in the library
whizard.a. Therefore, a main program written in Fortran 95 has
to contain the line
use whizard
in the preamble.25 For other
programming languages, the subroutines are accessed by a
platform-dependent prefix. For instance, the NAG compiler preprends
the string whizard_MP_ to all names in the module
whizard, so the subroutine whizard_integrate is
accessed by the name whizard_MP_whizard_integrate. The
access to arguments in subroutine calls may also be
platform-dependent.
The following routines are needed only if some non-standard behavior
of WHIZARD is intended. They closely match the command-line
arguments described in Sec. 5.6.
-
whizard_set_directory
- :
Set a working directory different from the default one.
-
Invocation:
-
call whizard_set_directory (directory)
- Arguments:
-
-
directory
- : string (input parameter)
- whizard_set_filename
- :
Set a specific filename to be used instead of
whizard for all
input and output files during this run. This is one way of avoiding
files getting overwritten. If an input file with specific filename does
not exist, the generic input files whizard.in etc. will be
used. The output files will always get the specific name, if any.
-
Invocation:
-
call whizard_set_filename (filename)
- Arguments:
-
-
filename
- : string (input parameter)
- whizard_set_input_file
- :
Define the name of the main input file. Needed only if it is
different from the default. The extension
.in will be appended.
-
Invocation:
-
call whizard_set_input_file (filename)
- Arguments:
-
-
filename
- : string (input parameter)
- whizard_enable_input_file
- :
Enable/disable reading of the main input file. If the logical
argument is
.false., calls to whizard_read_input
will not read any file. This makes sense if parameters are to be set by
different means, in particular by calls to
whizard_set_input.
-
Invocation:
-
call whizard_enable_input_file (flag)
- Arguments:
-
-
flag
- : logical true or false (input parameter)
- whizard_set_input
- :
Set input parameters directly. Each argument is a string of
parameter = value definitions in namelist format, as it would
appear in the input file enclosed by ``
&xxx'' and ``/'',
where xxx is the respective input block (e.g.,
process_input). The strings can be empty.
-
Invocation:
-
call whizard_set_input (process_input, integration_input,
simulation_input, diagnostics_input, parameter_input,
beam_input1, beam_input2)
- Arguments:
-
- process_input
- : namelist string (input parameter)
- integration_input
- : namelist string (input parameter)
- simulation_input
- : namelist string (input parameter)
- diagnostics_input
- : namelist string (input parameter)
- parameter_input
- : namelist string (input parameter)
- beam_input1
- : namelist string (input parameter)
- beam_input2
- : namelist string (input parameter)
The following routines are called, directly or indirectly, in standard
WHIZARD runs. The explicit interface allows for calling them
directly without reference to the main program.
-
whizard_read_input
- :
Start a new instance of WHIZARD and terminate any pending run.
Read the input file. If no specific filename has been given, this is
called
whizard.in. If this routine is not called explicitly,
it will be called automatically by whizard_integrate. Calling
it explicitly allows for overriding input data before they are fed
into the integration pass.
-
Invocation:
-
call whizard_read_input
- Arguments:
- none.
- whizard_number_processes
- :
Return the number of processes specified in the input file.
-
Invocation:
-
call whizard_number_processes (number)
- Arguments:
-
-
number
- : integer (output parameter)
- whizard_integrate
- : Calculate the total cross section
resp. read pre-generated results from file, as required by the input
data. This routine reads the model definition file
whizard.mdl, the cut configuration file
whizard.cut1 and, if required, it reads/writes the grids
whizard.grc and whizard.grb.
If this routine is not called explicitly, it will be called
automatically by the event generation routines.
-
Invocation:
-
call whizard_integrate
- Arguments:
- none.
- whizard_get_results
- : Read out the current
integration results for the indicated process: Total cross section,
error estimate, c2 value, reweighting efficiency estimate. If
process_index is zero, return the total integral, error and
efficiency average for all activated processes. The chi2 value
is zero in that case. If no integration has yet been performed, all
returned numbers are zero. If the process_index is larger than
the number of processes (as returned by
whizard_number_processes),
an error is issued.
The returned numbers are of the default real kind, which is double
precision under normal circumstances. The default real kind may be
accessed in a Fortran 95 program by a use kinds directive and
including the kinds-src subdirectory in the module search path.
-
Invocation:
-
call whizard_get_results (process_index, integral, error, chi2, efficiency)
- Arguments:
-
- process_index
- : integer (input variable)
- integral
- : default
real (output variable)
- error
- : default real (output
variable)
- chi2
- : default real (output variable)
- efficiency
- : default real (output variable)
- whizard_generate
- : Generate a fixed number of events
resp. read them from file, as specified by the input data. This
routine reads the cut/analysis configuration
whizard.cut5 and
reads/writes the event file(s) whizard.evt and
whizard.evx, if required. Furthermore, it can write histogram
data in whizard-plots.dat. This subroutine is equivalent
to consecutive calls of whizard_events_begin,
whizard_event (multiple times), and whizard_events_end
(see below).
-
Invocation:
-
call whizard_generate
- Arguments:
- none
- whizard_events_begin
- : Prepare event generation for a
n_calls events. It is assumed that the main
program contains a loop with exactly n_calls to
whizard_event. However, in case of unweighted events
generation, this number is used only for opening a STDHEP format file,
so the value is irrelevant if that option is not used.
-
Invocation:
-
call whizard_events_begin (n_calls)
- Arguments:
-
-
n_calls
- : 64-bit integer (input parameter)
- whizard_event
- : Generate a single (unweighted) event.
The event is stored in the standard
HEPEVT COMMON block. If
requested by the input data, it is also written to file. To make use of
the current event, the caller routine must be able to access
HEPEVT.
-
Invocation:
-
call whizard_event
- Arguments:
- none
- whizard_events_end
- : Finalize event generation.
-
Invocation:
-
call whizard_events_end
- Arguments:
- none
- whizard_end
- : Finalize the WHIZARD run. Close all
open files and deallocate memory.
-
Invocation:
-
call whizard_end
- Arguments:
- none