(git:3add494)
xas_tdp_atom Module Reference

This module deals with all the integrals done on local atomic grids in xas_tdp. This is mostly used to compute the xc kernel matrix elements wrt two RI basis elements (centered on the same excited atom) <P|fxc(r)|Q>, where the kernel fxc is purely a function of the ground state density and r. This is also used to compute the SOC matrix elements in the orbital basis. More...

Functions/Subroutines

subroutine, public init_xas_atom_env (xas_atom_env, xas_tdp_env, xas_tdp_control, qs_env, ltddfpt)
 Initializes a xas_atom_env type given the qs_enxas_atom_env, qs_envv. More...
 
subroutine, public init_xas_atom_grid_harmo (xas_atom_env, grid_info, do_xc, qs_env)
 Initializes the atomic grids and harmonics for the RI atomic calculations. More...
 
subroutine, public truncate_radial_grid (grid_atom, max_radius)
 Reduces the radial extension of an atomic grid such that it only covers a given radius. More...
 
subroutine, public compute_sphi_so (ikind, basis_type, sphi_so, qs_env)
 Computes the contraction coefficients to go from spherical orbitals to sgf for a given atomic kind and a given basis type. More...
 
subroutine, public calculate_density_coeffs (xas_atom_env, qs_env)
 Compute the coefficients to project the density on a partial RI_XAS basis. More...
 
subroutine, public integrate_fxc_atoms (int_fxc, xas_atom_env, xas_tdp_control, qs_env)
 Integrate the xc kernel as a function of r on the atomic grids for the RI_XAS basis. More...
 
subroutine, public integrate_soc_atoms (matrix_soc, xas_atom_env, qs_env, soc_atom_env)
 Computes the SOC matrix elements with respect to the ORB basis set for each atomic kind and put them as the block diagonal of dbcsr_matrix. More...
 

Detailed Description

This module deals with all the integrals done on local atomic grids in xas_tdp. This is mostly used to compute the xc kernel matrix elements wrt two RI basis elements (centered on the same excited atom) <P|fxc(r)|Q>, where the kernel fxc is purely a function of the ground state density and r. This is also used to compute the SOC matrix elements in the orbital basis.

Function/Subroutine Documentation

◆ init_xas_atom_env()

subroutine, public xas_tdp_atom::init_xas_atom_env ( type(xas_atom_env_type), pointer  xas_atom_env,
type(xas_tdp_env_type), pointer  xas_tdp_env,
type(xas_tdp_control_type), pointer  xas_tdp_control,
type(qs_environment_type), pointer  qs_env,
logical, optional  ltddfpt 
)

Initializes a xas_atom_env type given the qs_enxas_atom_env, qs_envv.

Parameters
xas_atom_envthe xas_atom_env to initialize
xas_tdp_env...
xas_tdp_control...
qs_env...
ltddfpt...

Definition at line 158 of file xas_tdp_atom.F.

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

◆ init_xas_atom_grid_harmo()

subroutine, public xas_tdp_atom::init_xas_atom_grid_harmo ( type(xas_atom_env_type), pointer  xas_atom_env,
character(len=default_string_length), dimension(:, :), pointer  grid_info,
logical, intent(in)  do_xc,
type(qs_environment_type), pointer  qs_env 
)

Initializes the atomic grids and harmonics for the RI atomic calculations.

Parameters
xas_atom_env...
grid_info...
do_xcWhether the xc kernel will ne computed on the atomic grids. If not, the harmonics are built for the orbital basis for all kinds.
qs_env...
Note
Largely inspired by init_rho_atom subroutine

Definition at line 236 of file xas_tdp_atom.F.

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

◆ truncate_radial_grid()

subroutine, public xas_tdp_atom::truncate_radial_grid ( type(grid_atom_type), pointer  grid_atom,
real(dp), intent(in)  max_radius 
)

Reduces the radial extension of an atomic grid such that it only covers a given radius.

Parameters
grid_atomthe atomic grid from which we truncate the radial part
max_radiusthe maximal radial extension of the resulting grid
Note
Since the RI density used for <P|fxc|Q> is only of quality close to the atom, and the integrand only non-zero within the radius of the gaussian P,Q. One can reduce the grid extansion to the largest radius of the RI basis set

Definition at line 402 of file xas_tdp_atom.F.

Here is the caller graph for this function:

◆ compute_sphi_so()

subroutine, public xas_tdp_atom::compute_sphi_so ( integer, intent(in)  ikind,
character(len=*), intent(in)  basis_type,
real(dp), dimension(:, :), pointer  sphi_so,
type(qs_environment_type), pointer  qs_env 
)

Computes the contraction coefficients to go from spherical orbitals to sgf for a given atomic kind and a given basis type.

Parameters
ikindthe kind for which we compute the coefficients
basis_typethe type of basis for which we compute
sphi_sowhere the new contraction coefficients are stored (not yet allocated)
qs_env...

Definition at line 449 of file xas_tdp_atom.F.

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

◆ calculate_density_coeffs()

subroutine, public xas_tdp_atom::calculate_density_coeffs ( type(xas_atom_env_type), pointer  xas_atom_env,
type(qs_environment_type), pointer  qs_env 
)

Compute the coefficients to project the density on a partial RI_XAS basis.

Parameters
xas_atom_env...
qs_env...
Note
The density is n = sum_ab P_ab*phi_a*phi_b, the RI basis covers the products of orbital sgfs => n = sum_ab sum_cd P_ab (phi_a phi_b xi_c) S_cd^-1 xi_d = sum_d coeff_d xi_d , where xi are the RI basis func. In this case, with the partial RI projection, the RI basis is restricted to an excited atom and its neighbors at a time. Leads to smaller overlap matrix to invert and less 3-center overlap to compute. The procedure is repeated for each excited atom

Definition at line 859 of file xas_tdp_atom.F.

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

◆ integrate_fxc_atoms()

subroutine, public xas_tdp_atom::integrate_fxc_atoms ( type(cp_2d_r_p_type), dimension(:, :), pointer  int_fxc,
type(xas_atom_env_type), pointer  xas_atom_env,
type(xas_tdp_control_type), pointer  xas_tdp_control,
type(qs_environment_type), pointer  qs_env 
)

Integrate the xc kernel as a function of r on the atomic grids for the RI_XAS basis.

Parameters
int_fxcthe global array containing the (P|fxc|Q) integrals, for all spin configurations
xas_atom_env...
xas_tdp_control...
qs_env...
Note
Note that if closed-shell, alpha-alpha term and beta-beta terms are the same Store the (P|fxc|Q) integrals on the processor they were computed on int_fxc(1)matrix is alpha-alpha, 2: alpha-beta, 3: beta-beta

Definition at line 1874 of file xas_tdp_atom.F.

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

◆ integrate_soc_atoms()

subroutine, public xas_tdp_atom::integrate_soc_atoms ( type(dbcsr_p_type), dimension(:), pointer  matrix_soc,
type(xas_atom_env_type), optional, pointer  xas_atom_env,
type(qs_environment_type), pointer  qs_env,
type(soc_atom_env_type), optional, pointer  soc_atom_env 
)

Computes the SOC matrix elements with respect to the ORB basis set for each atomic kind and put them as the block diagonal of dbcsr_matrix.

Parameters
matrix_socthe matrix where the SOC is stored
xas_atom_env...
qs_env...
soc_atom_env...
Note
We compute: <da_dx|V(4c^2-2V)|db_dy> - <da_dy|V(4c^2-2V)|db_dx>, where V is a model potential on the atomic grid. What we get is purely imaginary

Definition at line 3368 of file xas_tdp_atom.F.

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