ao_two_e_ints

Here, all two-electron integrals (\(1/r_{12}\)) are computed. As they have 4 indices and many are zero, they are stored in a map, as defined in utils/map_module.f90.

To fetch an AO integral, use the get_ao_two_e_integral(i,j,k,l,ao_integrals_map) function.

The conventions are: * For AO integrals : (ij|kl) = (11|22) = <ik|jl> = <12|12>

EZFIO parameters

io_ao_two_e_integrals

Read/Write AO integrals from/to disk [ Write | Read | None ]

Default: None

ao_integrals_threshold

If | (pq|rs) | < ao_integrals_threshold then (pq|rs) is zero

Default: 1.e-15

do_direct_integrals

Compute integrals on the fly (very slow, only for debugging)

Default: False

Providers

ao_integrals_cache

File : ao_two_e_ints/map_integrals.irp.f

double precision, allocatable   :: ao_integrals_cache   (0:64*64*64*64)

Cache of AO integrals for fast access

Needs:

  • ao_integrals_cache_min

  • ao_integrals_map

  • ao_two_e_integrals_in_map

ao_integrals_cache_max

File : ao_two_e_ints/map_integrals.irp.f

integer :: ao_integrals_cache_min
integer :: ao_integrals_cache_max

Min and max values of the AOs for which the integrals are in the cache

Needs:

  • ao_num

Needed by:

  • ao_integrals_cache

  • ao_integrals_cache_periodic

ao_integrals_cache_min

File : ao_two_e_ints/map_integrals.irp.f

integer :: ao_integrals_cache_min
integer :: ao_integrals_cache_max

Min and max values of the AOs for which the integrals are in the cache

Needs:

  • ao_num

Needed by:

  • ao_integrals_cache

  • ao_integrals_cache_periodic

ao_integrals_cache_periodic

File : ao_two_e_ints/map_integrals.irp.f

complex*16, allocatable :: ao_integrals_cache_periodic  (0:64*64*64*64)

Cache of AO integrals for fast access

Needs:

  • ao_integrals_cache_min

  • ao_integrals_map

  • ao_two_e_integrals_in_map

ao_integrals_map

File : ao_two_e_ints/map_integrals.irp.f

type(map_type)  :: ao_integrals_map

AO integrals

Needs:

  • ao_num

Needed by:

  • ao_integrals_cache

  • ao_integrals_cache_periodic

  • ao_two_e_integrals_in_map

  • mo_two_e_integral_jj_from_ao

  • mo_two_e_integrals_in_map

  • mo_two_e_integrals_vv_from_ao

ao_two_e_integral_schwartz

File : ao_two_e_ints/two_e_integrals.irp.f

double precision, allocatable   :: ao_two_e_integral_schwartz   (ao_num,ao_num)

Needed to compute Schwartz inequalities

Needs:

  • ao_coef_normalized_ordered_transp

  • ao_expo_ordered_transp

  • ao_nucl

  • ao_num

  • ao_power

  • ao_prim_num

  • n_pt_max_integrals

  • nucl_coord

ao_two_e_integrals_in_map

File : ao_two_e_ints/two_e_integrals.irp.f

logical :: ao_two_e_integrals_in_map
Map of Atomic integrals

i(r1) j(r2) 1/r12 k(r1) l(r2)

Needs:

  • ao_coef_normalized_ordered_transp

  • ao_expo_ordered_transp

  • ao_integrals_map

  • ao_nucl

  • ao_num

  • ao_power

  • ao_prim_num

  • ezfio_filename

  • io_ao_two_e_integrals

  • mpi_master

  • n_pt_max_integrals

  • nproc

  • nucl_coord

  • read_ao_two_e_integrals

  • zmq_context

  • zmq_socket_pull_tcp_address

  • zmq_state

Needed by:

  • ao_integrals_cache

  • ao_integrals_cache_periodic

  • mo_two_e_integral_jj_from_ao

  • mo_two_e_integrals_in_map

  • mo_two_e_integrals_vv_from_ao

gauleg_t2

File : ao_two_e_ints/gauss_legendre.irp.f

double precision, allocatable   :: gauleg_t2    (n_pt_max_integrals,n_pt_max_integrals/2)
double precision, allocatable   :: gauleg_w     (n_pt_max_integrals,n_pt_max_integrals/2)

t_w(i,1,k) = w(i) t_w(i,2,k) = t(i)

Needs:

  • n_pt_max_integrals

gauleg_w

File : ao_two_e_ints/gauss_legendre.irp.f

double precision, allocatable   :: gauleg_t2    (n_pt_max_integrals,n_pt_max_integrals/2)
double precision, allocatable   :: gauleg_w     (n_pt_max_integrals,n_pt_max_integrals/2)

t_w(i,1,k) = w(i) t_w(i,2,k) = t(i)

Needs:

  • n_pt_max_integrals

general_primitive_integral:

File : ao_two_e_ints/two_e_integrals.irp.f

  double precision function general_primitive_integral(dim,            &
P_new,P_center,fact_p,p,p_inv,iorder_p,                        &
Q_new,Q_center,fact_q,q,q_inv,iorder_q)

Computes the integral <pq|rs> where p,q,r,s are Gaussian primitives

Calls:

  • add_poly_multiply()

  • give_polynom_mult_center_x()

  • multiply_poly()

i_x1_new:

File : ao_two_e_ints/two_e_integrals.irp.f

recursive subroutine I_x1_new(a,c,B_10,B_01,B_00,res,n_pt)

recursive function involved in the two-electron integral

Needs:

  • n_pt_max_integrals

Called by:

  • i_x1_new()

  • i_x2_new()

  • integrale_new()

Calls:

  • i_x1_new()

  • i_x2_new()

i_x1_pol_mult_a1:

File : ao_two_e_ints/two_e_integrals.irp.f

recursive subroutine I_x1_pol_mult_a1(c,B_10,B_01,B_00,C_00,D_00,d,nd,n_pt_in)

Recursive function involved in the two-electron integral

Called by:

  • i_x1_pol_mult()

  • i_x1_pol_mult_a2()

  • i_x1_pol_mult_recurs()

Calls:

  • i_x2_pol_mult()

  • multiply_poly()

i_x1_pol_mult_a2:

File : ao_two_e_ints/two_e_integrals.irp.f

recursive subroutine I_x1_pol_mult_a2(c,B_10,B_01,B_00,C_00,D_00,d,nd,n_pt_in)

Recursive function involved in the two-electron integral

Called by:

  • i_x1_pol_mult()

  • i_x1_pol_mult_recurs()

Calls:

  • i_x1_pol_mult_a1()

  • i_x2_pol_mult()

  • multiply_poly()

i_x1_pol_mult_recurs:

File : ao_two_e_ints/two_e_integrals.irp.f

recursive subroutine I_x1_pol_mult_recurs(a,c,B_10,B_01,B_00,C_00,D_00,d,nd,n_pt_in)

Recursive function involved in the two-electron integral

Called by:

  • i_x1_pol_mult()

  • i_x1_pol_mult_recurs()

Calls:

  • i_x1_pol_mult_a1()

  • i_x1_pol_mult_a2()

  • i_x1_pol_mult_recurs()

  • multiply_poly()

i_x2_new:

File : ao_two_e_ints/two_e_integrals.irp.f

recursive subroutine I_x2_new(c,B_10,B_01,B_00,res,n_pt)

recursive function involved in the two-electron integral

Needs:

  • n_pt_max_integrals

Called by:

  • i_x1_new()

Calls:

  • i_x1_new()

i_x2_pol_mult:

File : ao_two_e_ints/two_e_integrals.irp.f

recursive subroutine I_x2_pol_mult(c,B_10,B_01,B_00,C_00,D_00,d,nd,dim)

Recursive function involved in the two-electron integral

Called by:

  • i_x1_pol_mult()

  • i_x1_pol_mult_a1()

  • i_x1_pol_mult_a2()

  • i_x2_pol_mult()

Calls:

  • i_x2_pol_mult()

  • multiply_poly()

Subroutines / functions

ao_idx2_sq:

File : ao_two_e_ints/map_integrals.irp.f

subroutine ao_idx2_sq(i,j,ij)

Called by:

  • two_e_integrals_index_2fold()

ao_idx2_sq_rev:

File : ao_two_e_ints/map_integrals.irp.f

subroutine ao_idx2_sq_rev(i,k,ik)

reverse square compound index

Called by:

  • two_e_integrals_index_reverse_2fold()

ao_idx2_tri_key:

File : ao_two_e_ints/map_integrals.irp.f

subroutine ao_idx2_tri_key(i,j,ij)

Called by:

  • two_e_integrals_index_2fold()

ao_idx2_tri_rev_key:

File : ao_two_e_ints/map_integrals.irp.f

subroutine ao_idx2_tri_rev_key(i,k,ik)

return i<=k

Called by:

  • two_e_integrals_index_reverse_2fold()

ao_l4:

File : ao_two_e_ints/two_e_integrals.irp.f

integer function ao_l4(i,j,k,l)

Computes the product of l values of i,j,k,and l

Needs:

  • ao_l

ao_two_e_integral:

File : ao_two_e_ints/two_e_integrals.irp.f

double precision function ao_two_e_integral(i,j,k,l)
integral of the AO basis <ik|jl> or (ij|kl)

i(r1) j(r1) 1/r12 k(r2) l(r2)

Needs:

  • ao_coef_normalized_ordered_transp

  • ao_expo_ordered_transp

  • ao_nucl

  • ao_power

  • ao_prim_num

  • n_pt_max_integrals

  • nucl_coord

Calls:

  • give_explicit_poly_and_gaussian()

ao_two_e_integral_schwartz_accel:

File : ao_two_e_ints/two_e_integrals.irp.f

double precision function ao_two_e_integral_schwartz_accel(i,j,k,l)
integral of the AO basis <ik|jl> or (ij|kl)

i(r1) j(r1) 1/r12 k(r2) l(r2)

Needs:

  • ao_coef_normalized_ordered_transp

  • ao_expo_ordered_transp

  • ao_integrals_threshold

  • ao_nucl

  • ao_power

  • ao_prim_num

  • n_pt_max_integrals

  • nucl_coord

Calls:

  • give_explicit_poly_and_gaussian()

ao_two_e_integral_zero:

File : ao_two_e_ints/screening.irp.f

logical function ao_two_e_integral_zero(i,j,k,l)

Needs:

  • ao_integrals_threshold

  • ao_overlap_abs

  • ao_two_e_integral_schwartz

  • is_periodic

  • read_ao_two_e_integrals

ao_two_e_integrals_in_map_collector:

File : ao_two_e_ints/integrals_in_map_slave.irp.f

subroutine ao_two_e_integrals_in_map_collector(zmq_socket_pull)

Collects results from the AO integral calculation

Needs:

  • ao_integrals_map

  • ao_num

Called by:

  • ao_two_e_integrals_in_map

Calls:

  • end_zmq_to_qp_run_socket()

  • insert_into_ao_integrals_map()

ao_two_e_integrals_in_map_slave:

File : ao_two_e_ints/integrals_in_map_slave.irp.f

subroutine ao_two_e_integrals_in_map_slave(thread,iproc)

Computes a buffer of integrals

Needs:

  • ao_num

Called by:

  • ao_two_e_integrals_in_map_slave_inproc()

  • ao_two_e_integrals_in_map_slave_tcp()

Calls:

  • compute_ao_integrals_jl()

  • end_zmq_push_socket()

  • end_zmq_to_qp_run_socket()

  • push_integrals()

ao_two_e_integrals_in_map_slave_inproc:

File : ao_two_e_ints/integrals_in_map_slave.irp.f

subroutine ao_two_e_integrals_in_map_slave_inproc(i)

Computes a buffer of integrals. i is the ID of the current thread.

Called by:

  • ao_two_e_integrals_in_map

Calls:

  • ao_two_e_integrals_in_map_slave()

ao_two_e_integrals_in_map_slave_tcp:

File : ao_two_e_ints/integrals_in_map_slave.irp.f

subroutine ao_two_e_integrals_in_map_slave_tcp(i)

Computes a buffer of integrals. i is the ID of the current thread.

Calls:

  • ao_two_e_integrals_in_map_slave()

clear_ao_map:

File : ao_two_e_ints/map_integrals.irp.f

subroutine clear_ao_map

Frees the memory of the AO map

Needs:

  • ao_integrals_map

Calls:

  • map_deinit()

compute_ao_integrals_jl:

File : ao_two_e_ints/two_e_integrals.irp.f

subroutine compute_ao_integrals_jl(j,l,n_integrals,buffer_i,buffer_value)

Parallel client for AO integrals

Needs:

  • ao_integrals_threshold

  • ao_num

Called by:

  • ao_two_e_integrals_in_map_slave()

Calls:

  • two_e_integrals_index()

compute_ao_two_e_integrals:

File : ao_two_e_ints/two_e_integrals.irp.f

subroutine compute_ao_two_e_integrals(j,k,l,sze,buffer_value)

Compute AO 1/r12 integrals for all i and fixed j,k,l

Needs:

  • ao_num

Called by:

  • mo_two_e_integral_jj_from_ao

  • mo_two_e_integrals_vv_from_ao

eri:

File : ao_two_e_ints/two_e_integrals.irp.f

double precision function ERI(alpha,beta,delta,gama,a_x,b_x,c_x,d_x,a_y,b_y,c_y,d_y,a_z,b_z,c_z,d_z)
ATOMIC PRIMTIVE two-electron integral between the 4 primitives ::

primitive_1 = x1**(a_x) y1**(a_y) z1**(a_z) exp(-alpha * r1**2) primitive_2 = x1**(b_x) y1**(b_y) z1**(b_z) exp(- beta * r1**2) primitive_3 = x2**(c_x) y2**(c_y) z2**(c_z) exp(-delta * r2**2) primitive_4 = x2**(d_x) y2**(d_y) z2**(d_z) exp(- gama * r2**2)

Calls:

  • integrale_new()

gauleg:

File : ao_two_e_ints/gauss_legendre.irp.f

subroutine gauleg(x1,x2,x,w,n)

Gauss-Legendre

Called by:

  • gauleg_t2

get_ao_map_size:

File : ao_two_e_ints/map_integrals.irp.f

function get_ao_map_size()

Returns the number of elements in the AO map

Needs:

  • ao_integrals_map

get_ao_two_e_integral:

File : ao_two_e_ints/map_integrals.irp.f

double precision function get_ao_two_e_integral(i,j,k,l,map) result(result)

Gets one AO bi-electronic integral from the AO map

Needs:

  • ao_integrals_cache

  • ao_integrals_cache_min

  • ao_two_e_integrals_in_map

Calls:

  • map_get()

  • two_e_integrals_index()

get_ao_two_e_integral_periodic:

File : ao_two_e_ints/map_integrals.irp.f

complex*16 function get_ao_two_e_integral_periodic(i,j,k,l,map) result(result)

Gets one AO bi-electronic integral from the AO map

Needs:

  • ao_integrals_cache_min

  • ao_integrals_cache_periodic

  • ao_integrals_map

  • ao_two_e_integrals_in_map

Calls:

  • map_get()

  • two_e_integrals_index_2fold()

get_ao_two_e_integrals:

File : ao_two_e_ints/map_integrals.irp.f

subroutine get_ao_two_e_integrals(j,k,l,sze,out_val)

Gets multiple AO bi-electronic integral from the AO map . All i are retrieved for j,k,l fixed. physicist convention : <ij|kl>

Needs:

  • ao_integrals_map

  • ao_two_e_integrals_in_map

Called by:

  • add_integrals_to_map()

  • add_integrals_to_map_no_exit_34()

  • add_integrals_to_map_three_indices()

get_ao_two_e_integrals_non_zero:

File : ao_two_e_ints/map_integrals.irp.f

subroutine get_ao_two_e_integrals_non_zero(j,k,l,sze,out_val,out_val_index,non_zero_int)

Gets multiple AO bi-electronic integral from the AO map . All non-zero i are retrieved for j,k,l fixed.

Needs:

  • ao_integrals_map

  • ao_integrals_threshold

  • ao_two_e_integrals_in_map

Called by:

  • mo_two_e_integral_jj_from_ao

  • mo_two_e_integrals_vv_from_ao

Calls:

  • map_get()

  • two_e_integrals_index()

get_ao_two_e_integrals_non_zero_jl:

File : ao_two_e_ints/map_integrals.irp.f

subroutine get_ao_two_e_integrals_non_zero_jl(j,l,thresh,sze_max,sze,out_val,out_val_index,non_zero_int)

Gets multiple AO bi-electronic integral from the AO map . All non-zero i are retrieved for j,k,l fixed.

Needs:

  • ao_integrals_map

  • ao_two_e_integrals_in_map

Calls:

  • map_get()

  • two_e_integrals_index()

get_ao_two_e_integrals_non_zero_jl_from_list:

File : ao_two_e_ints/map_integrals.irp.f

subroutine get_ao_two_e_integrals_non_zero_jl_from_list(j,l,thresh,list,n_list,sze_max,out_val,out_val_index,non_zero_int)

Gets multiple AO two-electron integrals from the AO map . All non-zero i are retrieved for j,k,l fixed.

Needs:

  • ao_integrals_map

  • ao_two_e_integrals_in_map

Calls:

  • map_get()

  • two_e_integrals_index()

get_ao_two_e_integrals_periodic:

File : ao_two_e_ints/map_integrals.irp.f

subroutine get_ao_two_e_integrals_periodic(j,k,l,sze,out_val)

Gets multiple AO bi-electronic integral from the AO map . All i are retrieved for j,k,l fixed. physicist convention : <ij|kl>

Needs:

  • ao_integrals_map

  • ao_two_e_integrals_in_map

give_polynom_mult_center_x:

File : ao_two_e_ints/two_e_integrals.irp.f

subroutine give_polynom_mult_center_x(P_center,Q_center,a_x,d_x,p,q,n_pt_in,pq_inv,pq_inv_2,p10_1,p01_1,p10_2,p01_2,d,n_pt_out)

subroutine that returns the explicit polynom in term of the “t” variable of the following polynomw :

$I_{x_1}(a_x,d_x,p,q) , I_{x_1}(a_y,d_y,p,q) I_{x_1}(a_z,d_z,p,q)$

Called by:

  • general_primitive_integral()

Calls:

  • i_x1_pol_mult()

i_x1_pol_mult:

File : ao_two_e_ints/two_e_integrals.irp.f

subroutine I_x1_pol_mult(a,c,B_10,B_01,B_00,C_00,D_00,d,nd,n_pt_in)

Recursive function involved in the two-electron integral

Called by:

  • give_polynom_mult_center_x()

Calls:

  • i_x1_pol_mult_a1()

  • i_x1_pol_mult_a2()

  • i_x1_pol_mult_recurs()

  • i_x2_pol_mult()

idx2_tri_int:

File : ao_two_e_ints/map_integrals.irp.f

subroutine idx2_tri_int(i,j,ij)
idx2_tri_rev_int:

File : ao_two_e_ints/map_integrals.irp.f

subroutine idx2_tri_rev_int(i,k,ik)

return i<=k

insert_into_ao_integrals_map:

File : ao_two_e_ints/map_integrals.irp.f

subroutine insert_into_ao_integrals_map(n_integrals,buffer_i, buffer_values)

Create new entry into AO map

Needs:

  • ao_integrals_map

Called by:

  • ao_two_e_integrals_in_map_collector()

Calls:

  • map_append()

integrale_new:

File : ao_two_e_ints/two_e_integrals.irp.f

subroutine integrale_new(I_f,a_x,b_x,c_x,d_x,a_y,b_y,c_y,d_y,a_z,b_z,c_z,d_z,p,q,n_pt)

Calculates the integral of the polynomial :

$I_{x_1}(a_x+b_x,c_x+d_x,p,q) , I_{x_1}(a_y+b_y,c_y+d_y,p,q) , I_{x_1}(a_z+b_z,c_z+d_z,p,q)$ in $( 0 ; 1)$

Needs:

  • gauleg_t2

  • n_pt_max_integrals

Called by:

  • eri()

Calls:

  • i_x1_new()

n_pt_sup:

File : ao_two_e_ints/two_e_integrals.irp.f

integer function n_pt_sup(a_x,b_x,c_x,d_x,a_y,b_y,c_y,d_y,a_z,b_z,c_z,d_z)

Returns the upper boundary of the degree of the polynomial involved in the two-electron integral :

$I_x(a_x,b_x,c_x,d_x) , I_y(a_y,b_y,c_y,d_y) , I_z(a_z,b_z,c_z,d_z)$

push_integrals:

File : ao_two_e_ints/integrals_in_map_slave.irp.f

subroutine push_integrals(zmq_socket_push, n_integrals, buffer_i, buffer_value, task_id)

Push integrals in the push socket

Called by:

  • ao_two_e_integrals_in_map_slave()

two_e_integrals_index:

File : ao_two_e_ints/map_integrals.irp.f

subroutine two_e_integrals_index(i,j,k,l,i1)

Gives a unique index for i,j,k,l using permtuation symmetry. i <-> k, j <-> l, and (i,k) <-> (j,l) for non-periodic systems

Called by:

  • ao_integrals_cache

  • ao_integrals_map

  • banned_excitation

  • compute_ao_integrals_jl()

  • four_idx_novvvv()

  • get_ao_two_e_integral()

  • get_ao_two_e_integrals_non_zero()

  • get_ao_two_e_integrals_non_zero_jl()

  • get_ao_two_e_integrals_non_zero_jl_from_list()

  • get_two_e_integral()

  • mo_integrals_cache

  • mo_integrals_map

two_e_integrals_index_2fold:

File : ao_two_e_ints/map_integrals.irp.f

subroutine two_e_integrals_index_2fold(i,j,k,l,i1)

Called by:

  • ao_integrals_cache_periodic

  • get_ao_two_e_integral_periodic()

Calls:

  • ao_idx2_sq()

  • ao_idx2_tri_key()

two_e_integrals_index_reverse:

File : ao_two_e_ints/map_integrals.irp.f

subroutine two_e_integrals_index_reverse(i,j,k,l,i1)

Computes the 4 indices $i,j,k,l$ from a unique index $i_1$. For 2 indices $i,j$ and $i le j$, we have $p = i(i-1)/2 + j$. The key point is that because $j < i$, $i(i-1)/2 < p le i(i+1)/2$. So $i$ can be found by solving $i^2 - i - 2p=0$. One obtains $i=1 + sqrt{1+8p}/2$ and $j = p - i(i-1)/2$. This rule is applied 3 times. First for the symmetry of the pairs (i,k) and (j,l), and then for the symmetry within each pair.

two_e_integrals_index_reverse_2fold:

File : ao_two_e_ints/map_integrals.irp.f

subroutine two_e_integrals_index_reverse_2fold(i,j,k,l,i1)

Calls:

  • ao_idx2_sq_rev()

  • ao_idx2_tri_rev_key()