(git:6a2e663)
xas_tdp_integrals Module Reference

3-center integrals machinery for the XAS_TDP method More...

Functions/Subroutines

subroutine, public create_pqx_tensor (pq_X, ab_nl, ac_nl, matrix_dist, blk_size_1, blk_size_2, blk_size_3, only_bc_same_center)
 Prepares a tensor to hold 3-center integrals (pq|X), where p,q are distributed according to the given 2d dbcsr distribution of the given . The third dimension of the tensor is iteslf not distributed (i.e. the t_pgrid's third dimension has size 1). The blocks are reserved according to the neighbor lists. More...
 
subroutine, public fill_pqx_tensor (pq_X, ab_nl, ac_nl, basis_set_list_a, basis_set_list_b, basis_set_list_c, potential_parameter, qs_env, only_bc_same_center, eps_screen)
 Fills the given 3 dimensional (pq|X) tensor with integrals. More...
 
subroutine, public build_xas_tdp_ovlp_nl (ab_list, basis_a, basis_b, qs_env, excited_atoms, ext_dist2d)
 Builds a neighbor lists set for overlaping 2-center S_ab, where b is restricted on a a given list of atoms. Used for Coulomb RI where (aI|P) = sum_b C_bI (ab|P), where contraction coeff C_bI is assumed to be non-zero only on excited atoms. More...
 
subroutine, public build_xas_tdp_3c_nl (ac_list, basis_a, basis_c, op_type, qs_env, excited_atoms, x_range, ext_dist2d)
 Builds a neighbor lists set taylored for 3-center integral within XAS TDP, such that only excited atoms are taken into account for the list_c. More...
 
subroutine, public get_opt_3c_dist2d (opt_3c_dist2d, ab_list, ac_list, basis_set_a, basis_set_b, basis_set_c, qs_env, only_bc_same_center)
 Returns an optimized distribution_2d for the given neighbor lists based on an evaluation of the cost of the corresponding 3-center integrals. More...
 
subroutine, public compute_ri_3c_exchange (ex_atoms, xas_tdp_env, xas_tdp_control, qs_env)
 Computes the RI exchange 3-center integrals (ab|c), where c is from the RI_XAS basis and centered on excited atoms and kind. The operator used is that of the RI metric. More...
 
subroutine, public compute_ri_3c_coulomb (xas_tdp_env, qs_env)
 Computes the RI Coulomb 3-center integrals (ab|c), where c is from the RI_XAS basis and centered on the excited atoms of xas_tdp_env. More...
 
subroutine, public compute_ri_exchange2_int (ex_kind, xas_tdp_env, xas_tdp_control, qs_env)
 Computes the two-center Exchange integral needed for the RI in kernel calculation. Stores the integrals in the xas_tdp_env as global (small) arrays. Does that for a given excited kind. The quantity stored is M^-1 (P|Q) M^-1, where M is the RI metric. If the metric is the same as the exchange potential, then we end up with the V-approximation (P|Q)^-1 By default (if no metric), the ri_m_potential is a copy of the x_potential. More...
 
subroutine, public compute_ri_coulomb2_int (ex_kind, xas_tdp_env, xas_tdp_control, qs_env)
 Computes the two-center Coulomb integral needed for the RI in kernel calculation. Stores the integrals (P|Q)^-1 in the xas_tdp_env as global (small) arrays. Does that for a given excited kind. More...
 

Detailed Description

3-center integrals machinery for the XAS_TDP method

Author
A. Bussy (03.2020)

Function/Subroutine Documentation

◆ create_pqx_tensor()

subroutine, public xas_tdp_integrals::create_pqx_tensor ( type(dbt_type), intent(out)  pq_X,
type(neighbor_list_set_p_type), dimension(:), pointer  ab_nl,
type(neighbor_list_set_p_type), dimension(:), pointer  ac_nl,
type(dbcsr_distribution_type), intent(in)  matrix_dist,
integer, dimension(:), intent(in)  blk_size_1,
integer, dimension(:), intent(in)  blk_size_2,
integer, dimension(:), intent(in)  blk_size_3,
logical, intent(in), optional  only_bc_same_center 
)

Prepares a tensor to hold 3-center integrals (pq|X), where p,q are distributed according to the given 2d dbcsr distribution of the given . The third dimension of the tensor is iteslf not distributed (i.e. the t_pgrid's third dimension has size 1). The blocks are reserved according to the neighbor lists.

Parameters
pq_Xthe tensor to store the integrals
ab_nlthe 1st and 2nd center neighbor list
ac_nlthe 1st and 3rd center neighbor list
matrix_dist...
blk_size_1the block size in the first dimension
blk_size_2the block size in the second dimension
blk_size_3the block size in the third dimension
only_bc_same_centeronly keep block if, for the corresponding integral (ab|c), b and c share the same center, i.e. r_bc = 0.0

Definition at line 126 of file xas_tdp_integrals.F.

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

◆ fill_pqx_tensor()

subroutine, public xas_tdp_integrals::fill_pqx_tensor ( type(dbt_type)  pq_X,
type(neighbor_list_set_p_type), dimension(:), pointer  ab_nl,
type(neighbor_list_set_p_type), dimension(:), pointer  ac_nl,
type(gto_basis_set_p_type), dimension(:), pointer  basis_set_list_a,
type(gto_basis_set_p_type), dimension(:), pointer  basis_set_list_b,
type(gto_basis_set_p_type), dimension(:), pointer  basis_set_list_c,
type(libint_potential_type)  potential_parameter,
type(qs_environment_type), pointer  qs_env,
logical, intent(in), optional  only_bc_same_center,
real(dp), intent(in), optional  eps_screen 
)

Fills the given 3 dimensional (pq|X) tensor with integrals.

Parameters
pq_Xthe tensor to fill
ab_nlthe neighbor list for the first 2 centers
ac_nlthe neighbor list for the first and third centers
basis_set_list_abasis sets for first center
basis_set_list_bbasis sets for second center
basis_set_list_cbasis sets for third center
potential_parameterthe operator for the integrals
qs_env...
only_bc_same_centersame as in create_pqX_tensor
eps_screenthreshold for possible screening
Note
The following indices are happily mixed within this routine: First center i,a,p Second center: j,b,q Third center: k,c,X

Definition at line 274 of file xas_tdp_integrals.F.

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

◆ build_xas_tdp_ovlp_nl()

subroutine, public xas_tdp_integrals::build_xas_tdp_ovlp_nl ( type(neighbor_list_set_p_type), dimension(:), pointer  ab_list,
type(gto_basis_set_p_type), dimension(:), pointer  basis_a,
type(gto_basis_set_p_type), dimension(:), pointer  basis_b,
type(qs_environment_type), pointer  qs_env,
integer, dimension(:), intent(in), optional  excited_atoms,
type(distribution_2d_type), optional, pointer  ext_dist2d 
)

Builds a neighbor lists set for overlaping 2-center S_ab, where b is restricted on a a given list of atoms. Used for Coulomb RI where (aI|P) = sum_b C_bI (ab|P), where contraction coeff C_bI is assumed to be non-zero only on excited atoms.

Parameters
ab_listthe neighbor list
basis_abasis set list for atom a
basis_bbasis set list for atom b
qs_env...
excited_atomsthe indices of the excited atoms on which b is centered
ext_dist2duse an external distribution2d

Definition at line 664 of file xas_tdp_integrals.F.

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

◆ build_xas_tdp_3c_nl()

subroutine, public xas_tdp_integrals::build_xas_tdp_3c_nl ( type(neighbor_list_set_p_type), dimension(:), pointer  ac_list,
type(gto_basis_set_p_type), dimension(:), pointer  basis_a,
type(gto_basis_set_p_type), dimension(:), pointer  basis_c,
integer, intent(in)  op_type,
type(qs_environment_type), pointer  qs_env,
integer, dimension(:), intent(in), optional  excited_atoms,
real(dp), intent(in), optional  x_range,
type(distribution_2d_type), optional, pointer  ext_dist2d 
)

Builds a neighbor lists set taylored for 3-center integral within XAS TDP, such that only excited atoms are taken into account for the list_c.

Parameters
ac_listthe neighbor list ready for 3-center integrals
basis_abasis set list for atom a
basis_cbasis set list for atom c
op_typeto indicate whther the list should be built with overlap, Coulomb or else in mind
qs_env...
excited_atomsthe indices of the excited atoms to consider (if not given, all atoms are taken)
x_rangein case some truncated/screened operator is used, gives its range
ext_dist2dexternal distribution_2d to be used
Note
Based on setup_neighbor_list with added features

Definition at line 763 of file xas_tdp_integrals.F.

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

◆ get_opt_3c_dist2d()

subroutine, public xas_tdp_integrals::get_opt_3c_dist2d ( type(distribution_2d_type), pointer  opt_3c_dist2d,
type(neighbor_list_set_p_type), dimension(:), pointer  ab_list,
type(neighbor_list_set_p_type), dimension(:), pointer  ac_list,
type(gto_basis_set_p_type), dimension(:), pointer  basis_set_a,
type(gto_basis_set_p_type), dimension(:), pointer  basis_set_b,
type(gto_basis_set_p_type), dimension(:), pointer  basis_set_c,
type(qs_environment_type), pointer  qs_env,
logical, intent(in), optional  only_bc_same_center 
)

Returns an optimized distribution_2d for the given neighbor lists based on an evaluation of the cost of the corresponding 3-center integrals.

Parameters
opt_3c_dist2dthe optimized distribution_2d
ab_list...
ac_list...
basis_set_a...
basis_set_b...
basis_set_c...
qs_env...
only_bc_same_center...

Definition at line 893 of file xas_tdp_integrals.F.

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

◆ compute_ri_3c_exchange()

subroutine, public xas_tdp_integrals::compute_ri_3c_exchange ( integer, dimension(:), intent(in)  ex_atoms,
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 RI exchange 3-center integrals (ab|c), where c is from the RI_XAS basis and centered on excited atoms and kind. The operator used is that of the RI metric.

Parameters
ex_atomsexcited atoms on which the third center is located
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
This routine is called once for each excited atom. Because there are many different a,b pairs involved, load balance is ok. This allows memory saving

Definition at line 1058 of file xas_tdp_integrals.F.

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

◆ compute_ri_3c_coulomb()

subroutine, public xas_tdp_integrals::compute_ri_3c_coulomb ( type(xas_tdp_env_type), pointer  xas_tdp_env,
type(qs_environment_type), pointer  qs_env 
)

Computes the RI Coulomb 3-center integrals (ab|c), where c is from the RI_XAS basis and centered on the excited atoms of xas_tdp_env.

Parameters
xas_tdp_env...
qs_env...
Note
The ri_3c_coul tensor of xas_tdp_env is defined and allocated here. Only computed once for the whole system (for optimized load balance). Ok because not too much memory needed

Definition at line 1145 of file xas_tdp_integrals.F.

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

◆ compute_ri_exchange2_int()

subroutine, public xas_tdp_integrals::compute_ri_exchange2_int ( integer, intent(in)  ex_kind,
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 two-center Exchange integral needed for the RI in kernel calculation. Stores the integrals in the xas_tdp_env as global (small) arrays. Does that for a given excited kind. The quantity stored is M^-1 (P|Q) M^-1, where M is the RI metric. If the metric is the same as the exchange potential, then we end up with the V-approximation (P|Q)^-1 By default (if no metric), the ri_m_potential is a copy of the x_potential.

Parameters
ex_kind...
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
Computes all these integrals in non-PBCs as we assume that the range is short enough that atoms do not exchange with their periodic images

Definition at line 1240 of file xas_tdp_integrals.F.

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

◆ compute_ri_coulomb2_int()

subroutine, public xas_tdp_integrals::compute_ri_coulomb2_int ( integer, intent(in)  ex_kind,
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 two-center Coulomb integral needed for the RI in kernel calculation. Stores the integrals (P|Q)^-1 in the xas_tdp_env as global (small) arrays. Does that for a given excited kind.

Parameters
ex_kind...
xas_tdp_env...
xas_tdp_control...
qs_env...

Definition at line 1399 of file xas_tdp_integrals.F.

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