(git:6a2e663)
xas_tdp_kernel Module Reference

All the kernel specific subroutines for XAS TDP calculations. More...

Functions/Subroutines

subroutine, public kernel_coulomb_xc (coul_ker, xc_ker, donor_state, xas_tdp_env, xas_tdp_control, qs_env)
 Computes, if asked for it, the Coulomb and XC kernel matrices, in the usuall matrix format. More...
 
subroutine, public kernel_exchange (ex_ker, donor_state, xas_tdp_env, xas_tdp_control, qs_env)
 Computes the exact exchange kernel matrix using RI. Returns an array of 2 matrices, which are: 1) the on-diagonal kernel: (ab|I_sigma J_tau) * delta_sigma,tau 2) the off-diagonal spin-conserving kernel: (aJ_sigma|I_tau b) * delta_sigma,tau An internal analysis determines which of the above are computed (can range from 0 to 2),. More...
 
subroutine, public reserve_contraction_blocks (matrices, ri_atom, qs_env)
 Reserves the blocks in of a dbcsr matrix as needed for RI 3-center contraction (aI|P) More...
 
subroutine, public contract2_ao_to_domo (contr_int, op_type, donor_state, xas_tdp_env, xas_tdp_control, qs_env)
 Contract the ri 3-center integrals stored in a tensor with repect to the donor MOs coeffs, for a given excited atom k => (aI|k) = sum_b c_Ib (ab|k) More...
 
subroutine, public ri_all_blocks_mm (contr_int, PQ)
 Multiply all the blocks of a contracted RI integral (aI|P) by a matrix of type (P|...|Q) More...
 

Detailed Description

All the kernel specific subroutines for XAS TDP calculations.

Author
A. Bussy (03.2019)

Function/Subroutine Documentation

◆ kernel_coulomb_xc()

subroutine, public xas_tdp_kernel::kernel_coulomb_xc ( type(dbcsr_type), intent(inout)  coul_ker,
type(dbcsr_p_type), dimension(:), pointer  xc_ker,
type(donor_state_type), pointer  donor_state,
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 
)

Computes, if asked for it, the Coulomb and XC kernel matrices, in the usuall matrix format.

Parameters
coul_kerpointer the the Coulomb kernel matrix (can be void pointer)
xc_kerarray of pointer to the different xc kernels (5 of them): 1) the restricted closed-shell singlet kernel 2) the restricted closed-shell triplet kernel 3) the spin-conserving open-shell xc kernel 4) the on-diagonal spin-flip open-shell xc kernel
donor_state...
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
Coulomb and xc kernel are put together in the same routine because they use the same RI Coulomb: (aI|Jb) = (aI|P) (P|Q)^-1 (Q|Jb) XC : (aI|fxc|Jb) = (aI|P) (P|Q)^-1 (Q|fxc|R) (R|S)^-1 (S|Jb) In the above formula, a,b label the sgfs The routine analyses the xas_tdp_control to know which kernel must be computed and how (open-shell, singlet, triplet, ROKS, LSD, etc...) On entry, the pointers should be allocated

Definition at line 80 of file xas_tdp_kernel.F.

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

◆ kernel_exchange()

subroutine, public xas_tdp_kernel::kernel_exchange ( type(dbcsr_p_type), dimension(:), pointer  ex_ker,
type(donor_state_type), pointer  donor_state,
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 
)

Computes the exact exchange kernel matrix using RI. Returns an array of 2 matrices, which are: 1) the on-diagonal kernel: (ab|I_sigma J_tau) * delta_sigma,tau 2) the off-diagonal spin-conserving kernel: (aJ_sigma|I_tau b) * delta_sigma,tau An internal analysis determines which of the above are computed (can range from 0 to 2),.

Parameters
ex_ker...
donor_state...
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
In the case of spin-conserving excitation, the kernel must later be multiplied by the usual Q projector. In the case of spin-flip, one needs to project the excitations coming from alpha donor MOs on the unoccupied beta MOs. This is done by multiplying by a Q projector where the alpha-alpha and beta-beta quadrants are swapped The ex_ker array should be allocated on entry (not the internals)

Definition at line 594 of file xas_tdp_kernel.F.

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

◆ reserve_contraction_blocks()

subroutine, public xas_tdp_kernel::reserve_contraction_blocks ( type(dbcsr_p_type), dimension(:), pointer  matrices,
integer, intent(in)  ri_atom,
type(qs_environment_type), pointer  qs_env 
)

Reserves the blocks in of a dbcsr matrix as needed for RI 3-center contraction (aI|P)

Parameters
matricesthe matrices for which blocks are reserved
ri_atomthe index of the atom on which RI is done (= all coeffs of I are there, and P too)
qs_env...
Note
the end product are normal type matrices that are possibly slightly spraser as matrix_s

Definition at line 902 of file xas_tdp_kernel.F.

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

◆ contract2_ao_to_domo()

subroutine, public xas_tdp_kernel::contract2_ao_to_domo ( type(dbcsr_p_type), dimension(:), pointer  contr_int,
character(len=*), intent(in)  op_type,
type(donor_state_type), pointer  donor_state,
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 
)

Contract the ri 3-center integrals stored in a tensor with repect to the donor MOs coeffs, for a given excited atom k => (aI|k) = sum_b c_Ib (ab|k)

Parameters
contr_intthe contracted integrals as array of dbcsr matrices
op_typefor which operator type we contract (COULOMB or EXCHANGE)
donor_state...
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
In the output matrices, (aI_b|k) is stored at block a,b where I_b is the partial contraction that only includes coeffs from atom b. Note that the contracted matrix is not symmetric. To get the fully contracted matrix over b, one need to add the block columns of (aI_b|k) (get an array of size nao*nsgfp). This step is unnessary in our case because we assume locality of donor state, and only one column od (aI_b|k) is pouplated

Definition at line 975 of file xas_tdp_kernel.F.

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

◆ ri_all_blocks_mm()

subroutine, public xas_tdp_kernel::ri_all_blocks_mm ( type(dbcsr_p_type), dimension(:)  contr_int,
real(dp), dimension(:, :), intent(in)  PQ 
)

Multiply all the blocks of a contracted RI integral (aI|P) by a matrix of type (P|...|Q)

Parameters
contr_intthe integral array
PQthe smaller matrix to multiply all blocks
Note
It is assumed that all non-zero blocks have the same number of columns. Can pass partial arrays, e.g. contr_int(1:3)

Definition at line 1266 of file xas_tdp_kernel.F.

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