Multidip
1.0
Multi-photon matrix elements
|
I/O routines used by MULTIDIP. More...
Data Types | |
type | kmatrix |
K-matrix file. More... | |
type | scatakcoeffs |
Photoionization wave function coefficients. More... | |
type | mappedmatrix |
Auxiliary data structure for matrix (potentially memory-mapped, or distributed) More... | |
type | moleculardata |
RMT molecular data file. More... | |
Functions/Subroutines | |
subroutine | load_mapped_matrix (this, filename, u, offset, rows, cols) |
Read or map a matrix from file. More... | |
subroutine | destruct_mapped_matrix (this) |
Finalize a (potentially mapped) matrix object. More... | |
subroutine | read_kmatrices (km, lukmt, nkset) |
Read K-matrix files. More... | |
subroutine | read_kmatrix_file (km) |
Read K-matrix file. More... | |
subroutine | reset_kmatrix_position (km, skip) |
Reset I/O pointer to start of K-matrices. More... | |
subroutine | get_kmatrix (km, kmat, skip) |
Read single K-matrix from the K-matrix file. More... | |
subroutine | read_wfncoeffs (ak, lusct) |
Read wave function coefficients from files. More... | |
subroutine | read_bndcoeffs (bnd, Ei, lubnd, mgvn0, stot0) |
Read inner bound wave function coefficients. More... | |
subroutine | read_wfncoeffs_file (ak) |
Read wave function coefficients for a single symmetry from a file. More... | |
subroutine | reset_wfncoeffs_position (ak, skip) |
Reset I/O pointer to start of Ak-coeffs. More... | |
subroutine | get_wfncoeffs (ak, E, Re_Ak, Im_Ak, skip) |
Read single Ak-matrix from the Ak-coeffs file. More... | |
subroutine | read_molecular_data (moldat, filename, mpiio, read_wmat2) |
Read RMT molecular_data file. 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, component, irrpair, transp, nf, nn, X, Y) |
Multiply by dipole matrix. More... | |
subroutine | assert_status (message, errcode) |
Check for error code. More... | |
integer function | clip (x, a, b) |
Clamp value in range. 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_energy_grid (escat) |
Write photoelectron energies to file. More... | |
subroutine | write_partial_wave_moments (moldat, M, nesc, suffix) |
Write partial wave moments. More... | |
subroutine | write_raw_dipoles (M, chains, nesc, stem) |
Write raw transition dipole moments. More... | |
subroutine | write_cross_section (cs, nesc, erange, filename) |
Write cross sections to a file. More... | |
Variables | |
integer | myproc = 1 |
integer | nprocs = 1 |
I/O routines used by MULTIDIP.
This module contains routines that read the necessary input files (K-matrices, scattering coefficients, molecular_data files) and return the needed subset of data.
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.
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 1126 of file multidip_io.F90.
subroutine multidip_io::apply_dipole_matrix | ( | type(moleculardata), intent(in) | moldat, |
integer, intent(in) | component, | ||
integer, intent(in) | irrpair, | ||
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.
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'.
The matrix argument X and Y have the following one-dimensional distribution among processes:
+---+---+---+---+ | | | | | | | | | | | 0 | 1 | 2 | 3 | | | | | | | | | | | +---+---+---+---+
However, when ScaLAPACK multiplication is desired, the arrays need to be redistributed into the standard two-dimensional block-cyclic form
+-+-+-+-+-+-+-+-+ |0|1|2|3|0|1|2|3| +-+-+-+-+-+-+-+-+ |1|2|3|0|1|2|3|0| +-+-+-+-+-+-+-+-+ |2|3|0|1|2|3|0|1| +-+-+-+-+-+-+-+-+ |3|0|1|2|3|0|1|2| +-+-+-+-+-+-+-+-+
and back. This is accomplished by the subroutine pdgemr2d
. Note however, that this subroutine is prone to failure when the matrices exceed 4-byte addressable number of elements. Working with such matrices then requires patched version of ScaLAPACK that use long integers internally.
Definition at line 984 of file multidip_io.F90.
subroutine multidip_io::assert_status | ( | character(len=*), intent(in) | message, |
integer, intent(in) | errcode | ||
) |
Check for error code.
If the provided error code is non-zero, print the error message and abort.
Definition at line 1090 of file multidip_io.F90.
integer function multidip_io::clip | ( | integer | x, |
integer | a, | ||
integer | b | ||
) |
Clamp value in range.
Definition at line 1108 of file multidip_io.F90.
subroutine multidip_io::destruct_mapped_matrix | ( | type(mappedmatrix), intent(inout) | this | ) |
Finalize a (potentially mapped) matrix object.
Safely deallocates the matrix. If it has been mapped, the pointer is only nullified; unmapping happens automatically in the desctructor of the file_mapping type.
Definition at line 277 of file multidip_io.F90.
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.
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 926 of file multidip_io.F90.
subroutine multidip_io::get_kmatrix | ( | class(kmatrix), intent(in) | km, |
real(wp), dimension(:, :), intent(inout), allocatable | kmat, | ||
integer, intent(in), optional | skip | ||
) |
Read single K-matrix from the K-matrix file.
Assuming that the K-matrix file associated with this object is correctly positioned, read the next K-matrix record. Symmetrize the K-matrix and store it into the allocatable two-dimensional array kmat
(re/allocate as necessary).
Definition at line 463 of file multidip_io.F90.
subroutine multidip_io::get_wfncoeffs | ( | class(scatakcoeffs), intent(in) | ak, |
real(wp), intent(inout) | E, | ||
real(wp), dimension(:, :), intent(inout), allocatable | Re_Ak, | ||
real(wp), dimension(:, :), intent(inout), allocatable | Im_Ak, | ||
integer, intent(in), optional | skip | ||
) |
Read single Ak-matrix from the Ak-coeffs file.
Assuming that the Ak-coeffs file associated with this object is correctly positioned, read the next Ak-matrix record. Also returns the scattering energy in atomic units.
Definition at line 682 of file multidip_io.F90.
subroutine multidip_io::load_mapped_matrix | ( | class(mappedmatrix), intent(inout) | this, |
character(len=*), intent(in) | filename, | ||
integer(int32), intent(in) | u, | ||
integer(int64), intent(in) | offset, | ||
integer(int32), intent(in) | rows, | ||
integer(int32), intent(in) | cols | ||
) |
Read or map a matrix from file.
Depending on the compilation settings, either map the given section of the file to memory, or simply allocate and read the chunk. Optionally, read and distribute the matrix in the block-cyclic way.
Definition at line 157 of file multidip_io.F90.
subroutine multidip_io::read_bndcoeffs | ( | real(wp), dimension(:), intent(out) | bnd, |
real(wp), intent(out) | Ei, | ||
integer, intent(in) | lubnd, | ||
integer, intent(in) | mgvn0, | ||
integer, intent(in) | stot0 | ||
) |
Read inner bound wave function coefficients.
Read the inner region expansion coefficients of the bound state as calculated by BOUND.
[out] | bnd | Vector of inner region expansion coefficients to fill. |
[out] | Ei | Bound state energy as calculated by BOUND. |
[in] | lubnd | File unit number with the (unformatted) BOUND output. |
[in] | mgvn0 | Expected MGVN (to check). |
[in] | stot0 | Expected total spin (to check). |
Definition at line 551 of file multidip_io.F90.
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.
Read all needed K-matrix files. These are needed to calculate photionization coefficients (Ak) during the calculation. Also perform consistency check between the energy samples in individual files.
Definition at line 307 of file multidip_io.F90.
subroutine multidip_io::read_kmatrix_file | ( | class(kmatrix), intent(inout) | km | ) |
Read K-matrix file.
Read metadata from a K-matrix file, count K-matrices and store the unit and dataset for later use in retrieval of the K-matrices themselves. The K-matrices are not being read into memory here to avoid exhausting RAM (particularly in parallel mode) if the K-matrices are large.
Definition at line 350 of file multidip_io.F90.
subroutine multidip_io::read_molecular_data | ( | type(moleculardata), intent(inout) | moldat, |
character(*), intent(in) | filename, | ||
logical, intent(in) | mpiio, | ||
logical, intent(in) | read_wmat2 | ||
) |
Read RMT molecular_data file.
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 726 of file multidip_io.F90.
subroutine multidip_io::read_wfncoeffs | ( | type(scatakcoeffs), dimension(:), intent(inout), allocatable | ak, |
integer, dimension(:), intent(in) | lusct | ||
) |
Read wave function coefficients from files.
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 525 of file multidip_io.F90.
subroutine multidip_io::read_wfncoeffs_file | ( | class(scatakcoeffs), intent(inout) | ak | ) |
Read wave function coefficients for a single symmetry from a file.
Definition at line 592 of file multidip_io.F90.
subroutine multidip_io::reset_kmatrix_position | ( | class(kmatrix), intent(in) | km, |
integer, intent(in), optional | skip | ||
) |
Reset I/O pointer to start of K-matrices.
Rewind the unit to the start of the associated dataset and read through to the very beginning of the actual K-matrix data. Also skip the given number of leading K-matrices. This prepares the file for the subsequent calls to get_kmatrix
.
Definition at line 428 of file multidip_io.F90.
subroutine multidip_io::reset_wfncoeffs_position | ( | class(scatakcoeffs), intent(in) | ak, |
integer, intent(in), optional | skip | ||
) |
Reset I/O pointer to start of Ak-coeffs.
Rewind the unit to the start of the associated dataset and read through to the very beginning of the actual Ak-coeffs data. Also skip the given number of leading Ak-coeffs. This prepares the file for the subsequent calls to get_wfncoeffs
.
Definition at line 647 of file multidip_io.F90.
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.
Multiply boundary amplitudes for each inner region state by an element of the provided vector 'v'. Writes the result into 'vw'.
Definition at line 1173 of file multidip_io.F90.
subroutine multidip_io::setup_angular_integrals | ( | type(moleculardata), intent(inout) | moldat | ) |
Store the angular integrals in a more convenient shape.
Convert the Gaunt integral storage from a list to an easily addressable matrix.
Definition at line 890 of file multidip_io.F90.
subroutine multidip_io::write_cross_section | ( | real(wp), dimension(:, :), intent(in) | cs, |
integer, intent(in) | nesc, | ||
integer, dimension(2), intent(in) | erange, | ||
character(len=*), intent(in) | filename | ||
) |
Write cross sections to a file.
Write cross sections to a file. The first column of the array is expected to be the photon energy.
Definition at line 1353 of file multidip_io.F90.
subroutine multidip_io::write_energy_grid | ( | real(wp), dimension(:), intent(in) | escat | ) |
Write photoelectron energies to file.
Write the column list of photoelectron scattering energies to file. Atomic units are used.
Definition at line 1204 of file multidip_io.F90.
subroutine multidip_io::write_partial_wave_moments | ( | type(moleculardata), intent(in) | moldat, |
complex(wp), dimension(:, :, :), intent(in), allocatable | M, | ||
integer, intent(in) | nesc, | ||
character(len=*), intent(in) | suffix | ||
) |
Write partial wave moments.
Produce tables with the multi-photon ionisation matrix elements per partial wave.
Definition at line 1223 of file multidip_io.F90.
subroutine multidip_io::write_raw_dipoles | ( | complex(wp), dimension(:, :, :, :), intent(in), allocatable | M, |
integer, dimension(:, :), intent(in), allocatable | chains, | ||
integer, intent(in) | nesc, | ||
character(len=*), intent(in) | stem | ||
) |
Write raw transition dipole moments.
Produce tables with the multi-photon ionisation matrix elements per partial wave and per ionisation chain.
Definition at line 1280 of file multidip_io.F90.
integer multidip_io::myproc = 1 |
Definition at line 48 of file multidip_io.F90.
integer multidip_io::nprocs = 1 |
Definition at line 49 of file multidip_io.F90.