Multidip  1.0
Multi-photon matrix elements
multidip_outer Module Reference

Routines related to outer region asymptotic channels. More...

Functions/Subroutines

subroutine test_outer_expansion (filename, moldat, irr, ck, ap, Etot)
 Evaluate wfn in channels for inside and outside of the R-matrix sphere. More...
 
subroutine evaluate_fundamental_solutions (moldat, r, irr, nopen, Ek, S, C, Sp, Cp, sqrtknorm)
 Evaluate asymptotic solutions at boundary. More...
 
subroutine calculate_k_matrix (moldat, nopen, irr, Etot, S, C, Sp, Cp, Kmat)
 Calculate (generalized) K-matrix. More...
 

Detailed Description

Routines related to outer region asymptotic channels.

Authors
J Benda
Date
2023

Function/Subroutine Documentation

◆ calculate_k_matrix()

subroutine multidip_outer::calculate_k_matrix ( type(moleculardata), intent(in)  moldat,
integer, intent(in)  nopen,
integer, intent(in)  irr,
real(wp), intent(in)  Etot,
real(wp), dimension(:), intent(in)  S,
real(wp), dimension(:), intent(in)  C,
real(wp), dimension(:), intent(in)  Sp,
real(wp), dimension(:), intent(in)  Cp,
real(wp), dimension(:, :), intent(inout)  Kmat 
)

Calculate (generalized) K-matrix.

Authors
J Benda
Date
2023

Alternative to RSOLVE. Calculates the K-matrix without propagation, simply by projecting the inner wave function on the asymptotic channels. Unlike RSOLVE, though, this subroutine can include the long-range dipole potential. Only permanent dipoles (diagonal couplings) are considered.

First, the R-matrix is calculated from boundary amplitudes and R-matrix poles. Next, the asymptotic standing wave solutions of the asymptotic potential (including monopole and dipole potentials) are evaluated using the subroutine COULCC. Asymptotic radial wave-function for a given partial wave is represented as a linear combination of Coulomb functions with generalized (generally complex) angular momenta that diagonalize the sum of the centrifugal and dipole potentials, neglecting the transition dipoles.

Finally, the generalized K-matrix is obtained by solution of the standard matrix equation

\[ (C - R C') K = R S' - S \,. \]

In this equation the diagonal matrices S and C consist of the regular and irregular solutions of the dipole-coupled outer region asymptotic equations.

Parameters
[in]moldatMolecular data class.
[in]nopenNumber of open channels.
[in]irrIrreducible representation index.
[in]EtotTotal energy of the whole system.
[in]SValue of the regular asymptotic solution per partial wave channel.
[in]CValue of the irregular asymptotic solution per partial wave channel.
[in]SpDerivative of the regular asymptotic solution per partial wave channel.
[in]CpDerivative of the irregular asymptotic solution per partial wave channel.
[in,out]KmatKmatrix to calculate.

Definition at line 199 of file multidip_outer.f90.

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

◆ evaluate_fundamental_solutions()

subroutine multidip_outer::evaluate_fundamental_solutions ( type(moleculardata), intent(in)  moldat,
real(wp), intent(in)  r,
integer, intent(in)  irr,
integer, intent(in)  nopen,
real(wp), dimension(:), intent(in)  Ek,
real(wp), dimension(:), intent(inout)  S,
real(wp), dimension(:), intent(inout)  C,
real(wp), dimension(:), intent(inout)  Sp,
real(wp), dimension(:), intent(inout)  Cp,
logical, intent(in), optional  sqrtknorm 
)

Evaluate asymptotic solutions at boundary.

Authors
J Benda
Date
2023

Obtain amplitudes and derivatives of the fundamental solutions at region boundary.

The evaluated functions use the same normalization as ASYWFN in cfasym.f. That is, the amplitudes contain an extra factor of 1/sqrt(k) in addition to the standard normalization of Coulomb functions. Similarly, the derivatives contain the same factor of 1/sqrt(k). Because the functions are differentiated with respect to 'r', this 1/sqrt(k) factor combines with the 'k' coming from derivative of the argument, giving sqrt(k) compared to standard normalization of the derivative of the Coulomb functions.

Parameters
[in]moldatMolecular data class.
[in]rEvaluation radius.
[in]irrIrreducible representation index.
[in]nopenRestrict evaluation of solutions to this number of (open or closed) partial wave channels.
[in]EkPhotoelectron kinetic energy in each channel.
[in,out]SRegular solution amplitude per partial wave channel.
[in,out]CIrregular solution amplitude per partial wave channel.
[in,out]SpRegular solution derivative per partial wave channel.
[in,out]CpIrregular solution derivative per partial wave channel.
[in]sqrtknormSet to .false. to avoid adding the 1/sqrt(k) factor discussed above.

Definition at line 124 of file multidip_outer.f90.

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

◆ test_outer_expansion()

subroutine multidip_outer::test_outer_expansion ( character(len=*), intent(in)  filename,
type(moleculardata), intent(in)  moldat,
integer, intent(in)  irr,
real(wp), dimension(:, :), intent(in)  ck,
real(wp), dimension(:, :), intent(in)  ap,
real(wp), intent(in)  Etot 
)

Evaluate wfn in channels for inside and outside of the R-matrix sphere.

Authors
J Benda
Date
2023

Obtain the single-electron radial wave function amplitudes in channels by projecting the inner wave function at all sampling radii contained in molecular_data. Continue these to twice the R-matrix radius by evaluating the outer expansion in the outer region. The outer expansion assumes only outgoing (or exponentially decaying real) radial functions.

The results are saved into a text file, which has the radial distance in the first column and sampled radial wave-functions in the remaining columns, one per partial wave channel. Every first column gives real part, every second the imaginary part.

Parameters
[in]filenameName of a text file to write.
[in]moldatMolecular data class.
[in]irrInde xof irreducible representation in molecular_data.
[in]ckInner region expansion coefficients (in terms of inner region eigenstates).
[in]apOuter region expansion coefficients (in terms of partial wave channel eigenfunctions).
[in]EtotTotal energy of the system.

Definition at line 54 of file multidip_outer.f90.

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