Previous Up Next
C Using O'Mega's Output
The structure of the outputfile, the calling convention and the required libraries depends on the target language, of course.
C.1 Fortran90/95
The Fortran95 module written by O'Mega has the following signature
module omega_amplitude
C.1.1 Libraries
The imported Fortran modules are
omega_kinds
defines omega_prec, which can be whatever the Fortran compiler supports. NB: the support libraries have not yet been tuned to give reliable answers for amplitudes with gauge cancellations in single precision.
omega95
defines the vertices for Dirac spinors in the chiral representation and vectors.
omega95_bispinors
is an alternative that defines the vertices for Dirac and Majorana spinors in the chiral representation and vectors using the Feynman rules of [14].
omega_parameters
defines the coupling constants
  use omega_kinds
  use omega95
  use omega_parameters
  implicit none
  private
C.1.2 Summary of Exported Functions
The functions and subroutines experted by the Fortran95 module are
C.1.3 Maintenance Functions
They currently do nothing, but are here for WHIZARD's [15] convenience
create
is called only once at the very beginning.
reset
is called whenever parameters are changed.
destroy
is called at most once at the very end.
  subroutine create ()
  end subroutine create
  subroutine reset ()
  end subroutine reset
  subroutine destroy ()
  end subroutine destroy
Allocate an array of the size used by the heuristic that suppresses vanishing helicity combinations
  interface allocate_zero
     module procedure allocate_zero_1, allocate_zero_2
  end interface
for join numbering of in and out states
  subroutine allocate_zero_1 (zero)
    integer, dimension(:,:), pointer :: zero
  end subroutine allocate_zero_index
and for separate numbering of in and out states
  subroutine allocate_zero_2 (zero)
    integer, dimension(:,:,:,:), pointer :: zero
  end subroutine allocate_zero_index_inout
C.1.4 Inquiry Functions
The total number of particles, the number of incoming particles and the number of outgoing particles:
  pure function number_particles () result (n)
    integer :: n
  end function number_particles
  pure function number_particles_in () result (n)
    integer :: n
  end function number_particles_in
  pure function number_particles_out () result (n)
    integer :: n
  end function number_particles_out
The spin states of all particles that can give non-zero results and their number. The tables are interpreted as
s(1:,i)
contains the helicities for each particle for the ith helicity combination.
  pure function number_spin_states () result (n)
    integer :: n
  end function number_spin_states
  pure subroutine spin_states (s)
    integer, dimension(:,:), intent(inout) :: s
  end subroutine spin_states
The spin states of the incoming particles that can give non-zero results and their number:
  pure function number_spin_states_in () result (n)
    integer :: n
  end function number_spin_states_in
  pure subroutine spin_states_in (s)
    integer, dimension(:,:), intent(inout) :: s
  end subroutine spin_states_in
The spin states of the outgoing particles that can give non-zero results and their number:
  pure function number_spin_states_out () result (n)
    integer :: n
  end function number_spin_states_out
  pure subroutine spin_states_out (s)
    integer, dimension(:,:), intent(inout) :: s
  end subroutine spin_states_out
The flavor combinations of all particles that can give non-zero results and their number. The tables are interpreted as
f(1:,i)
contains the PDG particle code for each particle for the ith helicity combination.
  pure function number_flavor_states () result (n)
    integer :: n
  end function number_flavor_states
  pure subroutine flavor_states (f)
    integer, dimension(:,:), intent(inout) :: f
  end subroutine flavor_states
The flavor combinations of the incoming particles that can give non-zero results and their number.
  pure function number_flavor_states_in () result (n)
    integer :: n
  end function number_flavor_states_in
  pure subroutine flavor_states_in (f)
    integer, dimension(:,:), intent(inout) :: f
  end subroutine flavor_states_in
The flavor combinations of the outgoing particles that can give non-zero results and their number.
  pure function number_flavor_states_out () result (n)
    integer :: n
  end function number_flavor_states_out
  pure subroutine flavor_states_out (f)
    integer, dimension(:,:), intent(inout) :: f
  end subroutine flavor_states_out
The flavor combinations of all particles that always can give a zero result and their number:
  pure function number_flavor_zeros () result (n)
    integer :: n
  end function number_flavor_zeros
  pure subroutine flavor_zeros (f)
    integer, dimension(:,:), intent(inout) :: f
  end subroutine flavor_zeros
The flavor combinations of the incoming particles that always can give a zero result and their number:
  pure function number_flavor_zeros_in () result (n)
    integer :: n
  end function number_flavor_zeros_in
  pure subroutine flavor_zeros_in (f)
    integer, dimension(:,:), intent(inout) :: f
  end subroutine flavor_zeros_in
The flavor combinations of the outgoing particles that always can give a zero result and their number:
  pure function number_flavor_zeros_out () result (n)
    integer :: n
  end function number_flavor_zeros_out
  pure subroutine flavor_zeros_out (f)
    integer, dimension(:,:), intent(inout) :: f
  end subroutine flavor_zeros_out
The same initial and final state can appear more than once in the tensor product and we must avoid double counting.
  pure subroutine multiplicities (a)
    integer, dimension(:), intent(inout) :: a
  end subroutine multiplicities
  pure subroutine multiplicities_in (a)
    integer, dimension(:), intent(inout) :: a
  end subroutine multiplicities_in
  pure subroutine multiplicities_out (a)
    integer, dimension(:), intent(inout) :: a
  end subroutine multiplicities_out

C.1.5 Amplitude
The function arguments of of the amplitude are
k(0:3,1:)
are the particle momenta: k(0:3,1) and k(0:3,2) are the incoming momenta, k(0:3,3:) are the outgoing momenta. All momenta are the physical momenta, i. e. forward time-like or light-like. The signs of the incoming momenta are flipped internally. Unless asked by a commandline parameter, O'Mega will not check the validity of the momenta.
s(1:)
are the helicities in the same order as the momenta. s=±1 signify s=±1/2 for fermions. s=0 makes no sense for fermions and massless vector bosons s=4 signifies an unphysical polarization for vector boson that the users are not supposed to use. Unless asked by a commandline parameter, O'Mega will not check the validity of the helicities.
f(1:)
are the PDG particle codes in the same order as the momenta.
  pure function amplitude (k, s, f) result (amp)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:), intent(in) :: s, f
    complex(kind=omega_prec) :: amp
  end function amplitude
Identical to amplitude (k, s, flavors(:,f)), where flavors has been filled by flavor_states:
  pure function amplitude_f (k, s, f) result (amp)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:), intent(in) :: s
    integer, intent(in) :: f
    complex(kind=omega_prec) :: amp
  end function amplitude_f
Identical to amplitude (k, spins(:,s), flavors(:,f)), where spins has been filled by spin_states and flavors has been filled by flavor_states:
  pure function amplitude_1 (k, s, f) result (amp)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, intent(in) :: s, f
    complex(kind=omega_prec) :: amp
  end function amplitude_1
Similar to amplitude_1, but with separate incoming and outgoing particles:
  pure function amplitude_2 &
       (k, s_in, f_in, s_out, f_out) result (amp)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, intent(in) :: s_in, f_in, s_out, f_out
    complex(kind=omega_prec) :: amp
  end function amplitude_2
The following are subroutines and not functions, since Fortran95 restricts arguments of pure functions to intent(in), but we need to update the counter for vanishing amplitudes.
zero(1:,1:)
an array containing the number of times a combination of spin index and flavor index yielded a vanishing amplitude. After a certain threshold, these combinations will be skipped. allocate_zero will allocate the correct size.
n
the current event count
  pure subroutine amplitude_nonzero (amp, k, s, f, zero, n)
    complex(kind=omega_prec), intent(out) :: amp
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:), intent(in) :: s, f
    integer, dimension(:,:), intent(inout) :: zero
    integer, intent(in) :: n
  end subroutine amplitude_nonzero
  pure subroutine amplitude_1_nonzero (amp, k, s, f, zero, n)
    complex(kind=omega_prec), intent(out) :: amp
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, intent(in) :: s, f
    integer, dimension(:,:), intent(inout) :: zero
    integer, intent(in) :: n
  end subroutine amplitude_1_nonzero
  pure subroutine amplitude_f_nonzero &
         (amp, k, s, f, zero, n)
    complex(kind=omega_prec), intent(out) :: amp
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:), intent(in) :: s
    integer, intent(in) :: f
    integer, dimension(:,:), intent(inout) :: zero
    integer, intent(in) :: n
  end subroutine amplitude_f_nonzero
zero(1:,1:,1:,1:)
an array containing the number of times a combination of incoming and outgoing spin indices and flavor indices yielded a vanishing amplitude. allocate_zero will allocate the correct size.
  pure subroutine amplitude_2_nonzero &
       (amp, k, s_in, f_in, s_out, f_out, zero, n)
    complex(kind=omega_prec), intent(out) :: amp
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, intent(in) :: s_in, f_in, s_out, f_out
    integer, dimension(:,:,:,:), intent(inout) :: zero
    integer, intent(in) :: n
  end subroutine amplitude_2_nonzero
  pure function symmetry (f) result (s)
    real(kind=omega_prec) :: s
    integer, dimension(:), intent(in) :: f
  end function symmetry

C.1.6 Summation
The the sums of squared matrix elements, the optional mask smask can be used to sum only a subset of helicities or flavors.
  pure function spin_sum_sqme (k, f, smask) result (amp2)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:), intent(in) :: f
    logical, dimension(:), intent(in), optional :: smask
    real(kind=omega_prec) :: amp2
  end function spin_sum_sqme
  pure function spin_sum_sqme_1 (k, f, smask) result (amp2)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, intent(in) :: f
    logical, dimension(:), intent(in), optional :: smask
    real(kind=omega_prec) :: amp2
  end function spin_sum_sqme_1
  pure function sum_sqme (k, smask, fmask) result (amp2)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    logical, dimension(:), intent(in), optional :: smask, fmask
    real(kind=omega_prec) :: amp2
  end function sum_sqme
  pure subroutine spin_sum_sqme_nonzero (amp2, k, f, zero, n, smask)
    real(kind=omega_prec), intent(out) :: amp2
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:), intent(in) :: f
    integer, dimension(:,:), intent(inout) :: zero
    integer, intent(in) :: n
    logical, dimension(:), intent(in), optional :: smask
  end subroutine spin_sum_sqme_nonzero
  pure subroutine spin_sum_sqme_1_nonzero (amp2, k, f, zero, n, smask)
    real(kind=omega_prec), intent(out) :: amp2
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, intent(in) :: f
    integer, dimension(:,:), intent(inout) :: zero
    integer, intent(in) :: n
    logical, dimension(:), intent(in), optional :: smask
  end subroutine spin_sum_sqme_1_nonzero
  pure subroutine sum_sqme_nonzero (amp2, k, zero, n, smask, fmask)
    real(kind=omega_prec), intent(out) :: amp2
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    integer, dimension(:,:), intent(inout) :: zero
    integer, intent(in) :: n
    logical, dimension(:), intent(in), optional :: smask, fmask
  end subroutine sum_sqme_masked_nonzero
C.1.7 Density Matrix Transforms
There are also utility functions that implement the transformation of density matrices directly
r ® r' = T r T     (18)
i. e.
r'ff' =
 
å
ii'
Tfi rii' Tf'i'*     (19)
and avoid double counting
  pure subroutine scatter_correlated (k, rho_in, rho_out)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    complex(kind=omega_prec), dimension(:,:,:,:), &
      intent(in) :: rho_in
    complex(kind=omega_prec), dimension(:,:,:,:), &
      intent(inout) :: rho_out
  end subroutine scatter_correlated
  pure subroutine scatter_correlated_nonzero &
         (k, rho_in, rho_out, zero, n)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    complex(kind=omega_prec), dimension(:,:,:,:), &
      intent(in) :: rho_in
    complex(kind=omega_prec), dimension(:,:,:,:), &
      intent(inout) :: rho_out
    integer, dimension(:,:,:,:), intent(inout) :: zero
    integer, intent(in) :: n
  end subroutine scatter_correlated_nonzero
In no off-diagonal density matrix elements of the initial state are known, the computation can be performed more efficiently:
r'f =
 
å
i
Tfi ri Tfi* =
 
å
i
|Tfi|2 ri     (20)
  pure subroutine scatter_diagonal (k, rho_in, rho_out)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    real(kind=omega_prec), dimension(:,:), intent(in) :: rho_in
    real(kind=omega_prec), dimension(:,:), intent(inout) :: rho_out
  end subroutine scatter_diagonal
  pure subroutine scatter_diagonal_nonzero &
         (k, rho_in, rho_out, zero, n)
    real(kind=omega_prec), dimension(0:,:), intent(in) :: k
    real(kind=omega_prec), dimension(:,:), intent(in) :: rho_in
    real(kind=omega_prec), dimension(:,:), intent(inout) :: rho_out
    integer, dimension(:,:,:,:), intent(inout) :: zero
    integer, intent(in) :: n
  end subroutine scatter_diagonal_nonzero
Finis.
end module omega_amplitude
NB: the name of the module can be changed by a commandline parameter and Fortran95 features like pure can be disabled as well.
C.2 FORTRAN77
The preparation of a FORTRAN77 target is straightforward, but tedious and will only be considered if there is sufficient demand and support.
C.3 HELAS
This target for the HELAS library [7] is incomplete and no longer maintained. It was used as an early benchmark for the Fortran90/95 library. No vector boson selfcouplings are supported.
C.4 C, C++ & Java
These targets does not exist yet and we solicit suggestions from C++ and Java experts on useful calling conventions and suppport libraries that blend well with the HEP environments based on these languages. At least one of the authors believes that Java would be a better choice, but the political momentum behind C++ might cause an early support for C++ anyway.
Previous Up Next