Multidip  1.0
Multi-photon matrix elements
multidip_io Module Reference

I/O routines used by MULTIDIP. More...

Data Types

type  KMatrix
 K-matrix file. More...
 
type  MolecularData
 RMT molecular data file. More...
 
type  ScatAkCoeffs
 Photoionization wave function coefficients. More...
 

Functions/Subroutines

subroutine read_kmatrices (km, lukmt, nkset)
 Read K-matrix files. More...
 
subroutine read_wfncoeffs (ak, lusct)
 Read wave function coefficients from file. More...
 
subroutine read_molecular_data (moldat)
 Read RMT molecular_data file. More...
 
subroutine destruct_molecular_data (moldat)
 Finalize the MolecularData object. More...
 
subroutine setup_angular_integrals (moldat)
 Store the angular integrals in a more convenient shape. More...
 
subroutine get_diptrans (moldat, I, iidip, ifdip)
 Return dipole transition descriptors. More...
 
subroutine apply_dipole_matrix (moldat, I, s, transp, nf, nn, X, Y)
 Multiply by dipole matrix. More...
 
subroutine apply_boundary_amplitudes (moldat, irr, transp, X, Y)
 Multiply by boundary amplitudes. More...
 
subroutine scale_boundary_amplitudes (moldat, irr, v, vw)
 Scale boundary amplitudes matrix by a diagonal matrix. More...
 
subroutine write_partial_wave_moments (moldat, M, suffix)
 Write partial wave moments. More...
 
subroutine write_cross_section (cs)
 Write cross section to a file. More...
 

Variables

integer myproc = 1
 
integer nprocs = 1
 
real(wp), parameter alpha = 1/137.03599907_wp
 Fine structure constant. More...
 
character(len=1), dimension(3), parameter compt = ['x', 'y', 'z']
 

Detailed Description

I/O routines used by MULTIDIP.

Author
J Benda
Date
2020

This module contains routines that read the necessary input files (K-matrices, scattering coefficients, molecular_data files) and return the needed subset of data.

Function/Subroutine Documentation

◆ apply_boundary_amplitudes()

subroutine multidip_io::apply_boundary_amplitudes ( type(moleculardata), intent(in)  moldat,
integer, intent(in)  irr,
character(len=1), intent(in)  transp,
real(wp), dimension(:,:), intent(inout)  X,
real(wp), dimension(:,:), intent(inout)  Y 
)

Multiply by boundary amplitudes.

Author
J Benda
Date
2020

Multiply real matrix X by the matrix of boundary amplitudes and store the result into the real matrix Y. When transp is 'T', then the elements of the matrix X are expected to correspond to outer region channels, while the elements of the matrix Y to the inner region states. Otherwise, the opposite is assumed.

Definition at line 542 of file multidip_io.F90.

Here is the caller graph for this function:

◆ apply_dipole_matrix()

subroutine multidip_io::apply_dipole_matrix ( type(moleculardata), intent(in)  moldat,
integer, intent(in)  I,
integer, intent(in)  s,
character(len=1), intent(in)  transp,
integer, intent(in)  nf,
integer, intent(in)  nn,
real(wp), dimension(:, :, :), intent(inout)  X,
real(wp), dimension(:, :, :), intent(inout)  Y 
)

Multiply by dipole matrix.

Author
J Benda
Date
2020

Multiply matrix 'X' by the dipole matrix corresponding to dipole component 'I' and transition index 's'. Store the result in matrix 'Y'. Use 'nn' rows of 'X', and 'nf' rows of 'Y'.

Definition at line 499 of file multidip_io.F90.

Here is the caller graph for this function:

◆ destruct_molecular_data()

subroutine multidip_io::destruct_molecular_data ( type(moleculardata), intent(inout)  moldat)

Finalize the MolecularData object.

Author
J Benda
Date
2020

This subroutine will release the file mapping objects (if used).

Definition at line 416 of file multidip_io.F90.

Here is the caller graph for this function:

◆ get_diptrans()

subroutine multidip_io::get_diptrans ( type(moleculardata), intent(in)  moldat,
integer, intent(in)  I,
integer, dimension(:), intent(inout), allocatable  iidip,
integer, dimension(:), intent(inout), allocatable  ifdip 
)

Return dipole transition descriptors.

Author
J Benda
Date
2020

Based on the Cartesian component index I = 1, 2, 3, (re)allocate the arrays 'iidip' and 'ifdip' and fill them with initial/final irreducible representations pairs for dipole matrices present in the molecular_data file structure 'moldat'.

Definition at line 468 of file multidip_io.F90.

Here is the caller graph for this function:

◆ read_kmatrices()

subroutine multidip_io::read_kmatrices ( type(kmatrix), dimension(:), intent(inout), allocatable  km,
integer, dimension(:), intent(in)  lukmt,
integer, dimension(:), intent(in)  nkset 
)

Read K-matrix files.

Author
J Benda
Date
2020

Read all needed K-matrix files. These are needed to calculate scattering coefficients (Ak) during the calculation.

Definition at line 128 of file multidip_io.F90.

Here is the caller graph for this function:

◆ read_molecular_data()

subroutine multidip_io::read_molecular_data ( type(moleculardata), intent(inout)  moldat)

Read RMT molecular_data file.

Author
J Benda
Date
2020

Read data from the RMT molecular_data file. This includes in particular both the inner and outer region transition dipole elements, as well as for instance Gaunt angular integrals for all partial waves used.

Definition at line 292 of file multidip_io.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_wfncoeffs()

subroutine multidip_io::read_wfncoeffs ( type(scatakcoeffs), dimension(:), intent(inout), allocatable  ak,
integer, dimension(:), intent(in)  lusct 
)

Read wave function coefficients from file.

Author
J Benda
Date
2020

Read wave function (Ak) coefficients from an unformatted file written by RSOLVE. This is an optional functionality meant for dabugging. It is recommended that the Ak coefficients are calculated by the present program instead to avoid storing large datasets on disk or in memory.

Definition at line 215 of file multidip_io.F90.

Here is the caller graph for this function:

◆ scale_boundary_amplitudes()

subroutine multidip_io::scale_boundary_amplitudes ( type(moleculardata), intent(in)  moldat,
integer, intent(in)  irr,
real(wp), dimension(:), intent(in), allocatable  v,
real(wp), dimension(:, :), intent(inout), allocatable  vw 
)

Scale boundary amplitudes matrix by a diagonal matrix.

Author
J Benda
Date
2020

Multiply boundary amplitudes for each inner region state by an element of the provided vector 'v'. Writes the result into 'vw'.

Definition at line 580 of file multidip_io.F90.

Here is the caller graph for this function:

◆ setup_angular_integrals()

subroutine multidip_io::setup_angular_integrals ( type(moleculardata), intent(inout)  moldat)

Store the angular integrals in a more convenient shape.

Author
J Benda
Date
2020

Convert the Gaunt integral storage from a list to an easily addressable matrix.

Definition at line 434 of file multidip_io.F90.

Here is the caller graph for this function:

◆ write_cross_section()

subroutine multidip_io::write_cross_section ( real(wp), dimension(:, :), intent(inout), allocatable  cs)

Write cross section to a file.

Author
J Benda
Date
2020

Write orientation-averaged cross section to a file.

Definition at line 652 of file multidip_io.F90.

Here is the caller graph for this function:

◆ write_partial_wave_moments()

subroutine multidip_io::write_partial_wave_moments ( type(moleculardata), intent(in)  moldat,
complex(wp), dimension(:, :, :), intent(inout), allocatable  M,
character(len=*), intent(in)  suffix 
)

Write partial wave moments.

Author
J Benda
Date
2020

Produce tables with the multi-photon ionisation matrix elements per partial wave.

Definition at line 605 of file multidip_io.F90.

Here is the caller graph for this function:

Variable Documentation

◆ alpha

real(wp), parameter multidip_io::alpha = 1/137.03599907_wp

Fine structure constant.

Definition at line 47 of file multidip_io.F90.

◆ compt

character(len=1), dimension(3), parameter multidip_io::compt = ['x', 'y', 'z']

Definition at line 49 of file multidip_io.F90.

◆ myproc

integer multidip_io::myproc = 1

Definition at line 44 of file multidip_io.F90.

◆ nprocs

integer multidip_io::nprocs = 1

Definition at line 45 of file multidip_io.F90.