Multidip  1.0
Multi-photon matrix elements
Multidip Documentation

Calculate multi-photon dipole elements

Author
J Benda
Date
2020

Description

Computes the multi-photon transition elements of a given order and the corresponding generalized cross section. Proceeeds in the following way (see also the call graph in multidip):

  1. The required data are read from specified K-matrix files (see multidip_io::read_kmatrices) and the RMT molecular_data file (see multidip_io::read_molecular_data). The K-matrix files should ideally contain the full K-matrices, i.e. including the closed channel contributions at the R-matrix radius. This can be achieved by setting IKTYPE = 1 in R-solve. When only the (default) open-open K-matrices are used, there is a possibility of below-threshold spurious spikes in the cross sections for multi-channel calculations.
  2. For all intermediate orders, the intermediate states are calculated at all scattering energies (see multidip_routines::solve_intermediate_state).
  3. The final stationary photoionization wave function is calculated at all scattering energies (see multidip_routines::extract_dipole_elements).
  4. Transition matrix elements and photoionisation cross section are obtained (see multidip_routines::calculate_observables).

While the implemented method works for an arbitrary perturbation order, there are some caveats that restrict the applicability of this program:

  1. Asymptotic series is used in place of Coulomb-Hankel functions, which does not work well close to thresholds.
  2. For (N + M)-photon ionisation, i.e. when M photons are absorbed in the continuum, the asymptotic complexity is proportional to O(M!). This practically restricts the tractable orders to at most 4 absorptions in the continuum.

Multidip can either reuse existing photoionization wave function coefficient files written by RSOLVE (see multidip_io::read_wfncoeffs) or calculate the wave function coefficients on the fly (see multidip_special::scatAk).

Command line options

The program understands the following command line options:

    --help          display brief help and exit
    --input FILE    use FILE as the input file (instead of standard input)
    --test          run the internal sanity checks and exit

Input parameters

When the command line option --input is not used, MULTIDIP will read the input parameters from the standard input. The input consists of the input namelist &mdip, which is expected to have the following form:

   &mdip
     order = 2,               ! two-photon transitions

     lusct = 800,801,802,803, ! photoionization Ak coefficient files to read (optional)
     lukmt = 190,191,192,194, ! K-matrix files to read (ideally, should be produced in RSOLVE with IKTYPE = 1)
     nkset = 1,1,1,1,         ! which set from the K-matrix file to use (default is to use the first set)

     polar(:,1) = (0,0),(0,0),(1,0),  ! polarisation of the first absorbed photon (default is 0, i.e. average over polarisations)
     polar(:,2) = (0,0),(0,0),(1,0),  ! polarisation of the second absorbed photon (default is 0, i.e. average over polarisations)

     omega(1) = 0,            ! fixed energy in a.u. of the first absorbed photon or 0 for auto
     omega(2) = 0,            ! fixed energy in a.u. of the second absorbed photon or 0 for auto

     verbose = .false.,       ! print detailed information about the execution
   /

The code will use the irreducible representations specified by the K-matrix files. These need to be present also in the file molecular_data. The user is responsible for providing all relevant K-matrix files (irreducible representations) for the specified photon polarisations. Otherwise, the results may be wrong.

The polarisation of \( k \)-th photon is given as three complex-valued components of a unit Cartesian vector: \( (\mathrm{Re}\,\epsilon_x^{(k)}, \mathrm{Im}\,\epsilon_x^{(k)}), (\mathrm{Re}\,\epsilon_y^{(k)}, \mathrm{Im}\,\epsilon_y^{(k)}), (\mathrm{Re}\,\epsilon_z^{(k)}, \mathrm{Im}\,\epsilon_z^{(k)}) \). When a polarisation vector has only real components (as in the above example), some compilers (like ifort, but not gfortran) will allow omission of the parentheses and specifying just the real part. When polarisation of some photon is specified as zero, no partial wave multipole moments will be written to files at the end of the calculation and the cross section will be averaged over polarisation of this photon.

Some fixed absorbed photon energies can be specified on input using the omega parameter. There must always be at least one photon with variable energy (corresponds to omega = 0). The energies of the variable-energy photons will be then equally assigned so that the total energy of the system corresponds to the total energies on input (in the K-matrix files). By default omega = 0 for all photons, so all photons have variable energy and will correspond to the same wave length.

Output files

In all cases, the program produces the file gen_photo_xsec.txt with the generalized cross section. The first column is the energy of the first photon in eV, the second column is the total cross section in atomic units (area^N time^(N-1)). See multidip_routines::calculate_observables for details. The cross section is averaged over polarisations that are not given in the input file.

If all photon polarisations are given, the full partial wave N-photon transition matrix elements will be written to files pwdip-M-(I,L,M).txt and pwdip-M-(I,L,M)+cphase.txt, where M is the irreducible representation (0-based), I is the (1-based) index of the residual ion, and L and M are the angular momentum quantum numbers of the emitted electron's partial wave. The file with "+cphase" ending contains the additional partial wave phase factors i^L exp(-i sigma). The first column of these files is the real part of the matrix element, the second is the imaginary part. The total energies are not present, but correspond to the same energies as given by the input K-matrix files.

Configuration options

CMake option default value
WITH_COARRAYS OFF
WITH_GSL OFF
WITH_MMAP OFF

This program can optionally use special functions (the complex gamma function, Coulomb functions and the confluent hypergeometric functions) from the GNU Scientific Library instead of the UKRmol-out library. For that, pass the option -D WITH_GSL=ON to CMake during configuration.

This program is not explicitly parallelized, i.e. there are no direct calls to MPI and no OpenMP regions. However, it does benefit from optimized and/or threaded BLAS. Moreover, when compiled with the coarray support (and the flag -D WITH_COARRAYS=ON is passed to CMake), the loop over energies will be distributed among the parallel images in the round-robin fashion.

When configured with the -D WITH_MMAP=ON option, the inner region dipole matrices will not be read into memory, but only mapped into virtual memory. (This requires presence of such functionality in GBTOlib.)