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_pw_transition_elements).

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_routines::apply_Ak_coefficients_multidip).

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

     rmt_data = 'molecular_data',     ! name of the molecular RMT data file (default is 'molecular_data')

     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)
     lubnd = 500,                     ! unit with inner region expansion coefficients produced by BOUND (optional)

     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 (default is 0)
     omega(2) = 0,                    ! fixed energy in a.u. of the second absorbed photon or 0 for auto (default is 0)

     first_IP = 0,                    ! override ground state energy to yield given IP in a.u. (default is 0 = not used)
     verbose = .false.,               ! print detailed information about the execution (default is .false.)
     mpiio = .false.,                 ! run in HPC mode (requires MPI and ScaLAPACK, default is .false.)
   /

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.

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

When the polarisations of the photons are given, the program produces the file gen_photo_xsec.txt with the generalized cross section for oriented molecule. 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)), further columns are partial cross sections associated with transition to a specific final ionic state. See multidip_routines::calculate_pw_transition_elements for details.

Likewise, 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 \( \mathrm{i}^\ell \exp(-\mathrm{i} \sigma_\ell) \). 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.

Irrespective of the polarisations in the input namelist, MULTIDIP will produce the files gen_<n>photo_beta_<L>.txt with the parameters \( \beta_L \) of the molecular-orientatio-averaged photoelectron angular distribution,

\[ \frac{\mathrm{d}\sigma^{(p_1\dots p_n)}}{\mathrm{d}\Omega} = \frac{1}{4\pi} \beta_0 \left(1 + \beta_1 P_1(\cos\theta) + \dots + \beta_{2n} P_{2n}(\cos\theta) \right) , \]

where \( \beta_0 = \sigma^{(p_1\dots p_n)} \) is the integrated partial cross section and \( \beta_{L > 0} \) is the dimensionless asymmetry parameter for angular momentum \( L \). Currently, the laboratory-frame polarisations \( p_1,\dots,p_n \) are all hardcoded to zero (linear polarisation). Files that would contain zeros only are not written. The layout of these files is the same as of gen_photo_xsec.txt.

Finally, if some fully custom processing of the raw transition dipoles is required, it is possible to let them MULTIDIP write all to disk by use of the keyword raw. When set to the string "xyz" or "sph" or "both", raw transition dipoles in the selected angular basis will be written to files rawdip-CHAIN-(I,L,M).txt. Here "CHAIN" will be the absorption history of the transition, for instance "xx" if that raw dipole corresponds to absorption of two x-polarized photons in the Cartesian basis, or "mp0" if it corresponds to absorption of \( m = 0 \), \( m = +1 \) and \( m = -1 \) photons in the spherical basis. The first absorbed photon corresponds to the right-most index, as is customary in the typical bra-ket notation. The files contain real and imaginary parts of the dipoles for all scattering energies. The utility programs multidip_smooth_rawdips and multidip_rabitt_delays work with these files.

Configuration options

CMake option default value
WITH_GSL OFF
WITH_MMAP OFF
SCALAPACK_LIBRARIES (empty)

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.

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.)

This program benefits from optimized and/or threaded BLAS. Some sections of the code are also explicitly threaded (using OpenMP). Moreover, when compiled with the MPI support (this is controlled by GBTOlib level), the loop over energies (as well as some other pieces of the code) will be distributed among the parallel images in the round-robin fashion.

When the support for parallel MPI is enabled during compilation, it is also possible to enable use of MPI-IO and ScaLAPACK. MPI-IO is used to read the large inner region dipole matrices. ScaLAPACK is used to redistribute solutions and to multiply them with the inner dipoles. This MPI-IO/ScaLAPACK support will be compiled in if the option SCALAPACK_LIBRARIES is passed to CMake, and the code will be activated by setting mpiio = .true. in the program's input namelist.