utils

Contains general purpose utilities (sorting, maps, etc).

EZFIO parameters

restore_symm

If true, try to find symmetry in the MO coefficient matrices

Default: False

Providers

au_to_d

File : utils/units.irp.f

double precision        :: ha_to_ev
double precision        :: au_to_d
double precision        :: planck_cte
double precision        :: light_speed
double precision        :: ha_to_j
double precision        :: ha_to_nm

Some conversion between different units

binom

File : utils/util.irp.f

double precision, allocatable   :: binom        (0:40,0:40)
double precision, allocatable   :: binom_transp (0:40,0:40)

Binomial coefficients

Needed by:

  • binom_int

  • dettocsftransformationmatrix

  • nsomomax

  • psi_csf_coef

binom_int

File : utils/util.irp.f

integer*8, allocatable  :: binom_int    (0:40,0:40)
integer*8, allocatable  :: binom_int_transp     (0:40,0:40)

Binomial coefficients, as integers*8

Needs:

  • binom

Needed by:

  • dominant_dets_of_cfgs

  • n_dominant_dets_of_cfgs

  • psi_configuration_to_psi_det

binom_int_transp

File : utils/util.irp.f

integer*8, allocatable  :: binom_int    (0:40,0:40)
integer*8, allocatable  :: binom_int_transp     (0:40,0:40)

Binomial coefficients, as integers*8

Needs:

  • binom

Needed by:

  • dominant_dets_of_cfgs

  • n_dominant_dets_of_cfgs

  • psi_configuration_to_psi_det

binom_transp

File : utils/util.irp.f

double precision, allocatable   :: binom        (0:40,0:40)
double precision, allocatable   :: binom_transp (0:40,0:40)

Binomial coefficients

Needed by:

  • binom_int

  • dettocsftransformationmatrix

  • nsomomax

  • psi_csf_coef

degree_max_integration_lebedev

File : utils/angular_integration.irp.f

integer :: degree_max_integration_lebedev

integrate correctly a polynom of order “degree_max_integration_lebedev” needed for the angular integration according to LEBEDEV formulae

Needed by:

  • n_points_integration_angular_lebedev

  • theta_angular_integration_lebedev

dtranspose:

File : utils/transpose.irp.f

recursive subroutine dtranspose(A,LDA,B,LDB,d1,d2)

Transpose input matrix A into output matrix B

Called by:

  • dtranspose()

  • h_s2_u_0_nstates_openmp()

  • h_s2_u_0_nstates_zmq()

  • h_s2_u_0_two_e_nstates_openmp()

  • h_u_0_nstates_openmp()

  • h_u_0_nstates_zmq()

  • orb_range_2_rdm_openmp()

  • orb_range_2_rdm_state_av_openmp()

  • orb_range_2_trans_rdm_openmp()

Calls:

  • dtranspose()

fact_inv

File : utils/util.irp.f

double precision, allocatable   :: fact_inv     (128)

1/n!

give_explicit_cpoly_and_cgaussian:

File : utils/cgtos_utils.irp.f

subroutine give_explicit_cpoly_and_cgaussian(P_new, P_center, p, fact_k, iorder, &
       alpha, beta, a, b, Ae_center, Be_center, Ap_center, Bp_center, dim)

Transforms the product of

(x - x_Ap)^a(1) (x - x_Bp)^b(1) exp(-alpha (x - x_Ae)^2) exp(-beta (x - x_Be)^2) x (y - y_Ap)^a(2) (y - y_Bp)^b(2) exp(-alpha (y - y_Ae)^2) exp(-beta (y - y_Be)^2) x (z - z_Ap)^a(3) (z - z_Bp)^b(3) exp(-alpha (z - z_Ae)^2) exp(-beta (z - z_Be)^2)

into
fact_k * [sum (l_x = 0,i_order(1)) P_new(l_x,1) * (x-P_center(1))^l_x] exp (-p (x-P_center(1))^2)
  • [sum (l_y = 0,i_order(2)) P_new(l_y,2) * (y-P_center(2))^l_y] exp (-p (y-P_center(2))^2)

  • [sum (l_z = 0,i_order(3)) P_new(l_z,3) * (z-P_center(3))^l_z] exp (-p (z-P_center(3))^2)

WARNING ::: IF fact_k is too smal then: returns a “s” function centered in zero with an inifinite exponent and a zero polynom coef

Called by:

  • ao_2e_cgtos_schwartz_accel()

  • ao_two_e_integral_cgtos()

  • overlap_cgaussian_xyz()

Calls:

  • cgaussian_product()

  • multiply_cpoly()

  • recentered_cpoly2()

give_explicit_cpoly_and_cgaussian_x:

File : utils/cgtos_utils.irp.f

subroutine give_explicit_cpoly_and_cgaussian_x(P_new, P_center, p, fact_k, iorder, &
                                       alpha, beta, a, b, Ae_center, Be_center, Ap_center, Bp_center, dim)

Transform the product of

(x - x_Ap)^a (x - x_Bp)^b exp(-alpha (r - Ae)^2) exp(-beta (r - Be)^2)

into

fact_k sum_{i=0}^{iorder} (x - x_P)^i exp(-p (r - P)^2)

Called by:

  • overlap_cgaussian_x()

Calls:

  • multiply_cpoly()

  • recentered_cpoly2()

ha_to_ev

File : utils/units.irp.f

double precision        :: ha_to_ev
double precision        :: au_to_d
double precision        :: planck_cte
double precision        :: light_speed
double precision        :: ha_to_j
double precision        :: ha_to_nm

Some conversion between different units

ha_to_j

File : utils/units.irp.f

double precision        :: ha_to_ev
double precision        :: au_to_d
double precision        :: planck_cte
double precision        :: light_speed
double precision        :: ha_to_j
double precision        :: ha_to_nm

Some conversion between different units

ha_to_nm

File : utils/units.irp.f

double precision        :: ha_to_ev
double precision        :: au_to_d
double precision        :: planck_cte
double precision        :: light_speed
double precision        :: ha_to_j
double precision        :: ha_to_nm

Some conversion between different units

inv_int

File : utils/util.irp.f

double precision, allocatable   :: inv_int      (128)

1/i

light_speed

File : utils/units.irp.f

double precision        :: ha_to_ev
double precision        :: au_to_d
double precision        :: planck_cte
double precision        :: light_speed
double precision        :: ha_to_j
double precision        :: ha_to_nm

Some conversion between different units

n_points_integration_angular_lebedev

File : utils/angular_integration.irp.f

integer :: n_points_integration_angular_lebedev

Number of points needed for the angular integral

Needs:

  • degree_max_integration_lebedev

Needed by:

  • theta_angular_integration_lebedev

nproc

File : utils/util.irp.f

integer :: nproc

Number of current OpenMP threads

Needed by:

  • ao_two_e_integrals_erf_in_map

  • ao_two_e_integrals_in_map

  • cholesky_ao_num

  • h_apply_buffer_allocated

  • n_det

  • nthreads_davidson

  • nthreads_pt2

overlap_cgaussian_xyz:

File : utils/cgtos_one_e.irp.f

subroutine overlap_cgaussian_xyz(Ae_center, Be_center, alpha, beta, power_A, power_B, &
                         Ap_center, Bp_center, overlap_x, overlap_y, overlap_z, overlap, dim)
S_x = int (x - Ap_x)^{a_x} exp(-alpha (x - Ae_x)^2)

(x - Bp_x)^{b_x} exp(-beta (x - Be_x)^2) dx

S = S_x S_y S_z

for complex arguments

Called by:

  • ao_coef_norm_cgtos

  • ao_deriv2_cgtos_x

  • ao_overlap_cgtos

Calls:

  • give_explicit_cpoly_and_cgaussian()

phi_angular_integration_lebedev

File : utils/angular_integration.irp.f

double precision, allocatable   :: theta_angular_integration_lebedev    (n_points_integration_angular_lebedev)
double precision, allocatable   :: phi_angular_integration_lebedev      (n_points_integration_angular_lebedev)
double precision, allocatable   :: weights_angular_integration_lebedev  (n_points_integration_angular_lebedev)

Theta phi values together with the weights values for the angular integration : integral [dphi,dtheta] f(x,y,z) = 4 * pi * sum (1<i<n_points_integration_angular_lebedev) f(xi,yi,zi) Note that theta and phi are in DEGREES !!

Needs:

  • degree_max_integration_lebedev

  • n_points_integration_angular_lebedev

planck_cte

File : utils/units.irp.f

double precision        :: ha_to_ev
double precision        :: au_to_d
double precision        :: planck_cte
double precision        :: light_speed
double precision        :: ha_to_j
double precision        :: ha_to_nm

Some conversion between different units

qp_max_mem

File : utils/memory.irp.f

integer :: qp_max_mem

Maximum memory in Gb

Needs:

  • file_lock

  • mpi_master

Needed by:

  • ao_two_e_integral_alpha_chol

  • cholesky_ao_num

  • mo_two_e_integrals_erf_in_map

  • mo_two_e_integrals_in_map

  • pt2_j

  • pt2_w

shiftfact_op5_inv

File : utils/util.irp.f

double precision, allocatable   :: shiftfact_op5_inv    (128)

1 / Gamma(n + 0.5)

theta_angular_integration_lebedev

File : utils/angular_integration.irp.f

double precision, allocatable   :: theta_angular_integration_lebedev    (n_points_integration_angular_lebedev)
double precision, allocatable   :: phi_angular_integration_lebedev      (n_points_integration_angular_lebedev)
double precision, allocatable   :: weights_angular_integration_lebedev  (n_points_integration_angular_lebedev)

Theta phi values together with the weights values for the angular integration : integral [dphi,dtheta] f(x,y,z) = 4 * pi * sum (1<i<n_points_integration_angular_lebedev) f(xi,yi,zi) Note that theta and phi are in DEGREES !!

Needs:

  • degree_max_integration_lebedev

  • n_points_integration_angular_lebedev

transpose:

File : utils/transpose.irp.f

recursive subroutine transpose(A,LDA,B,LDB,d1,d2)

Transpose input matrix A into output matrix B

Called by:

  • transpose()

Calls:

  • transpose()

weights_angular_integration_lebedev

File : utils/angular_integration.irp.f

double precision, allocatable   :: theta_angular_integration_lebedev    (n_points_integration_angular_lebedev)
double precision, allocatable   :: phi_angular_integration_lebedev      (n_points_integration_angular_lebedev)
double precision, allocatable   :: weights_angular_integration_lebedev  (n_points_integration_angular_lebedev)

Theta phi values together with the weights values for the angular integration : integral [dphi,dtheta] f(x,y,z) = 4 * pi * sum (1<i<n_points_integration_angular_lebedev) f(xi,yi,zi) Note that theta and phi are in DEGREES !!

Needs:

  • degree_max_integration_lebedev

  • n_points_integration_angular_lebedev

Subroutines / functions

add_cpoly:

File : utils/cgtos_utils.irp.f

subroutine add_cpoly(b, nb, c, nc, d, nd)

Add two complex polynomials D(t) =! D(t) +( B(t) + C(t))

add_cpoly_multiply:

File : utils/cgtos_utils.irp.f

subroutine add_cpoly_multiply(b, nb, cst, d, nd)

Add a complex polynomial multiplied by a complex constant D(t) =! D(t) +( cst * B(t))

Called by:

  • general_primitive_integral_cgtos()

add_poly:

File : utils/integration.irp.f

subroutine add_poly(b,nb,c,nc,d,nd)

Add two polynomials D(t) =! D(t) +( B(t)+C(t))

add_poly_multiply:

File : utils/integration.irp.f

subroutine add_poly_multiply(b,nb,cst,d,nd)

Add a polynomial multiplied by a constant D(t) =! D(t) +( cst * B(t))

Called by:

  • general_primitive_integral()

  • general_primitive_integral_erf()

apply_rotation:

File : utils/linear_algebra.irp.f

subroutine apply_rotation(A,LDA,R,LDR,B,LDB,m,n)

Apply the rotation found by find_rotation

Calls:

  • dgemm()

approx_dble:

File : utils/util.irp.f

double precision function approx_dble(a,n)
binom_func:

File : utils/util.irp.f

double precision function binom_func(i,j)
cgaussian_product:

File : utils/cgtos_utils.irp.f

subroutine cgaussian_product(a, xa, b, xb, k, p, xp)

complex Gaussian product e^{-a (r-r_A)^2} e^{-b (r-r_B)^2} = k e^{-p (r-r_P)^2}

Called by:

  • give_explicit_cpoly_and_cgaussian()

cgaussian_product_x:

File : utils/cgtos_utils.irp.f

subroutine cgaussian_product_x(a, xa, b, xb, k, p, xp)

complex Gaussian product in 1D. e^{-a (x-x_A)^2} e^{-b (x-x_B)^2} = K e^{-p (x-x_P)^2}

check_mem:

File : utils/memory.irp.f

subroutine check_mem(rss_in,routine)

Checks if n gigabytes can be allocated. If not, exit the run.

Needs:

  • qp_max_mem

Called by:

  • ao_two_e_integral_alpha_chol

  • create_selection_buffer()

  • dav_double_dressed()

  • davidson_diag_csf_hjj()

  • davidson_diag_hjj()

  • davidson_diag_hjj_sjj()

  • davidson_diag_nonsym_hjj()

  • davidson_general()

  • davidson_general_diag_dressed_ext_rout_nonsym_b1space()

  • davidson_general_ext_rout()

  • davidson_general_ext_rout_diag_dressed()

  • davidson_general_ext_rout_dressed()

  • davidson_general_ext_rout_nonsym_b1space()

  • make_selection_buffer_s2()

  • pt2_collector()

  • pt2_j

  • pt2_w

  • remove_duplicates_in_selection_buffer()

  • run_cipsi()

  • run_slave_main()

  • run_stochastic_cipsi()

  • selection_collector()

  • testteethbuilding()

  • zmq_pt2()

Calls:

  • print_memory_usage()

  • resident_memory()

check_sym:

File : utils/util.irp.f

subroutine check_sym(A, n)
cpx_erf:

File : utils/cpx_erf.irp.f

complex*16 function cpx_erf(x, y)

compute erf(z) for z = x + i y

REF: Abramowitz and Stegun

cpx_erf_1:

File : utils/cpx_erf.irp.f

complex*16 function cpx_erf_1(x, y)

compute erf(z) for z = x + i y

REF: Abramowitz and Stegun

crint:

File : utils/cgtos_utils.irp.f

complex*16 function crint(n, rho)
crint_1:

File : utils/cpx_boys.irp.f

complex*16 function crint_1(n, rho)

Calls:

  • zboysfun00_1()

crint_1_vec:

File : utils/cpx_boys.irp.f

subroutine crint_1_vec(n_max, rho, vals)

Called by:

  • crint_sum()

Calls:

  • crint_smallz_vec()

  • zboysfun00_1()

crint_2:

File : utils/cpx_boys.irp.f

complex*16 function crint_2(n, rho)

Calls:

  • zboysfun()

  • zboysfunnrp()

crint_2_vec:

File : utils/cpx_boys.irp.f

subroutine crint_2_vec(n_max, rho, vals)

Calls:

  • crint_smallz_vec()

  • zboysfun()

  • zboysfunnrp()

crint_quad_1:

File : utils/cpx_boys.irp.f

subroutine crint_quad_1(n, rho, n_quad, crint_quad)
crint_quad_12:

File : utils/cpx_boys.irp.f

subroutine crint_quad_12(n, rho, n_quad, crint_quad)

Called by:

  • crint_quad_12_vec()

crint_quad_12_vec:

File : utils/cpx_boys.irp.f

subroutine crint_quad_12_vec(n_max, rho, vals)

Calls:

  • crint_quad_12()

crint_quad_2:

File : utils/cpx_boys.irp.f

subroutine crint_quad_2(n, rho, n_quad, crint_quad)
crint_smallz:

File : utils/cpx_boys.irp.f

complex*16 function crint_smallz(n, rho)

Standard version of rint

crint_smallz_vec:

File : utils/cpx_boys.irp.f

subroutine crint_smallz_vec(n_max, rho, vals)

Standard version of rint

Called by:

  • crint_1_vec()

  • crint_2_vec()

crint_sum:

File : utils/cpx_boys.irp.f

complex*16 function crint_sum(n_pt_out, rho, d1)

Calls:

  • crint_1_vec()

dble_fact:

File : utils/util.irp.f

double precision function dble_fact(n)
dble_fact_even:

File : utils/util.irp.f

double precision function dble_fact_even(n) result(fact2)

n!!

dble_fact_odd:

File : utils/util.irp.f

double precision function dble_fact_odd(n) result(fact2)

n!!

dble_logfact:

File : utils/util.irp.f

double precision function dble_logfact(n) result(logfact2)

n!!

derf_mu_x:

File : utils/util.irp.f

double precision function derf_mu_x(mu,x)
diag_mat_per_fock_degen:

File : utils/block_diag_degen.irp.f

subroutine diag_mat_per_fock_degen(fock_diag, mat_ref, n, thr_d, thr_nd, thr_deg, leigvec, reigvec, eigval)

subroutine that diagonalizes a matrix mat_ref BY BLOCK

the blocks are defined by the elements having the SAME DEGENERACIES in the entries “fock_diag”

examples : all elements having degeneracy 1 in fock_diag (i.e. not being degenerated) will be treated together

: all elements having degeneracy 2 in fock_diag (i.e. two elements are equal) will be treated together

: all elements having degeneracy 3 in fock_diag (i.e. two elements are equal) will be treated together

etc… the advantage is to guarentee no spurious mixing because of numerical problems.

Calls:

  • dsort()

  • give_degen_full_list()

  • isort()

  • non_hrmt_bieig()

diag_mat_per_fock_degen_core:

File : utils/block_diag_degen_core.irp.f

subroutine diag_mat_per_fock_degen_core(fock_diag, mat_ref, listcore,ncore, n, thr_d, thr_nd, thr_deg, leigvec, reigvec, eigval)

subroutine that diagonalizes a matrix mat_ref BY BLOCK

the blocks are defined by the elements having the SAME DEGENERACIES in the entries “fock_diag”

the elements of listcore are untouched

examples : all elements having degeneracy 1 in fock_diag (i.e. not being degenerated) will be treated together

: all elements having degeneracy 2 in fock_diag (i.e. two elements are equal) will be treated together

: all elements having degeneracy 3 in fock_diag (i.e. two elements are equal) will be treated together

etc… the advantage is to guarentee no spurious mixing because of numerical problems.

Calls:

  • dsort()

  • give_degen_full_listcore()

  • isort()

  • non_hrmt_bieig()

diag_nonsym_right:

File : utils/linear_algebra.irp.f

subroutine diag_nonsym_right(n, A, A_ldim, V, V_ldim, energy, E_ldim)

Called by:

  • davidson_diag_nonsym_hjj()

  • davidson_general_diag_dressed_ext_rout_nonsym_b1space()

  • davidson_general_ext_rout_nonsym_b1space()

Calls:

  • dgeevx()

  • dsort()

diagonalize_sym_matrix:

File : utils/util.irp.f

subroutine diagonalize_sym_matrix(N, A, e)

Diagonalize a symmetric matrix

Calls:

  • dsyev()

dset_order:

File : utils/sort.irp.f_template_90

subroutine dset_order(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A.

Called by:

  • ao_coef_norm_cgtos_ord

  • ao_coef_normalized_ordered

  • h_s2_u_0_nstates_openmp()

  • h_s2_u_0_nstates_zmq()

  • h_s2_u_0_two_e_nstates_openmp()

  • h_u_0_nstates_openmp()

  • h_u_0_nstates_zmq()

  • orb_range_2_rdm_openmp()

  • orb_range_2_rdm_state_av_openmp()

  • orb_range_2_trans_rdm_openmp()

  • psi_bilinear_matrix_transp_values

  • psi_bilinear_matrix_values

  • restore_symmetry()

dset_order_big:

File : utils/sort.irp.f_template_90

subroutine dset_order_big(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A. This is a version for very large arrays where the indices need to be in integer*8 format

eigsvd:

File : utils/linear_algebra.irp.f

subroutine eigSVD(A,LDA,U,LDU,D,Vt,LDVt,m,n)

Algorithm 3 of https://arxiv.org/pdf/1810.06860.pdf

A(m,n) = U(m,n) D(n) Vt(n,n) with m>n

Calls:

  • dgemm()

  • dscal()

  • lapack_diagd()

  • svd()

erf_e:

File : utils/cpx_erf.irp.f

complex*16 function erf_E(x, yabs)
erf_f:

File : utils/cpx_erf.irp.f

double precision function erf_F(x, yabs)
erf_g:

File : utils/cpx_erf.irp.f

complex*16 function erf_G(x, yabs)
erf_h:

File : utils/cpx_erf.irp.f

complex*16 function erf_H(x, yabs)
exp_matrix:

File : utils/linear_algebra.irp.f

subroutine exp_matrix(X,n,exp_X)

exponential of the matrix X: X has to be ANTI HERMITIAN !!

taken from Hellgaker, jorgensen, Olsen book

section evaluation of matrix exponential (Eqs. 3.1.29 to 3.1.31)

Calls:

  • dgemm()

  • get_a_squared()

  • lapack_diagd()

exp_matrix_taylor:

File : utils/linear_algebra.irp.f

subroutine exp_matrix_taylor(X,n,exp_X,converged)

exponential of a general real matrix X using the Taylor expansion of exp(X)

returns the logical converged which checks the convergence exponential of X using Taylor expansion

Called by:

  • umat

Calls:

  • dgemm()

extrapolate_data:

File : utils/extrapolation.irp.f

subroutine extrapolate_data(N_data, data, pt2, output)

Extrapolate the data to the FCI limit

Called by:

  • increment_n_iter()

Calls:

  • get_pseudo_inverse()

f_integral:

File : utils/integration.irp.f

double precision function F_integral(n,p)

function that calculates the following integral int_{-infty}^{+infty} x^n exp(-p x^2) dx

fact:

File : utils/util.irp.f

double precision function fact(n)

n!

fc_integral:

File : utils/cgtos_utils.irp.f

complex*16 function Fc_integral(n, inv_sq_p)

function that calculates the following integral int_{-infty}^{+infty} x^n exp(-p x^2) dx for complex valued p

find_rotation:

File : utils/linear_algebra.irp.f

subroutine find_rotation(A,LDA,B,m,C,n)

Find A.C = B

Calls:

  • dgemm()

  • get_pseudo_inverse()

format_w_error:

File : utils/format_w_error.irp.f

subroutine format_w_error(value,error,size_nb,max_nb_digits,format_value,str_error)

Format for double precision, value(error)

Called by:

  • pt2_collector()

Calls:

  • lock_io()

  • unlock_io()

gaussian_product:

File : utils/integration.irp.f

subroutine gaussian_product(a,xa,b,xb,k,p,xp)

Gaussian product in 1D. e^{-a (x-x_A)^2} e^{-b (x-x_B)^2} = K_{ab}^x e^{-p (x-x_P)^2}

Called by:

  • give_explicit_poly_and_gaussian()

  • give_explicit_poly_and_gaussian_double()

gaussian_product_v:

File : utils/integration.irp.f

subroutine gaussian_product_v(a, xa, LD_xa, b, xb, k, p, xp, n_points)

Gaussian product in 1D. e^{-a (x-x_A)^2} e^{-b (x-x_B)^2} = K_{ab}^x e^{-p (x-x_P)^2}

Using multiple A centers

Called by:

  • give_explicit_poly_and_gaussian_v()

gaussian_product_x:

File : utils/integration.irp.f

subroutine gaussian_product_x(a,xa,b,xb,k,p,xp)

Gaussian product in 1D. e^{-a (x-x_A)^2} e^{-b (x-x_B)^2} = K_{ab}^x e^{-p (x-x_P)^2}

Called by:

  • overlap_gaussian_xyz()

gaussian_product_x_v:

File : utils/integration.irp.f

subroutine gaussian_product_x_v(a,xa,b,xb,k,p,xp,n_points)

Gaussian product in 1D with multiple xa e^{-a (x-x_A)^2} e^{-b (x-x_B)^2} = K_{ab}^x e^{-p (x-x_P)^2}

get_a_squared:

File : utils/linear_algebra.irp.f

subroutine get_A_squared(A,n,A2)

A2 = A A where A is n x n matrix. Use the dgemm routine

Called by:

  • exp_matrix()

Calls:

  • dgemm()

get_ab_prod:

File : utils/linear_algebra.irp.f

subroutine get_AB_prod(A,n,m,B,l,AB)

AB = A B where A is n x m, B is m x l. Use the dgemm routine

Calls:

  • dgemm()

get_inverse:

File : utils/linear_algebra.irp.f

subroutine get_inverse(A,LDA,m,C,LDC)

Returns the inverse of the square matrix A

Called by:

  • ao_ortho_canonical_coef_inv

  • overlap_states

Calls:

  • dgetrf()

  • dgetri()

get_inverse_complex:

File : utils/linear_algebra.irp.f

subroutine get_inverse_complex(A,LDA,m,C,LDC)

Returns the inverse of the square matrix A

Calls:

  • zgetrf()

  • zgetri()

get_pseudo_inverse:

File : utils/linear_algebra.irp.f

subroutine get_pseudo_inverse(A, LDA, m, n, C, LDC, cutoff)

Find C = A^-1

Called by:

  • ao_cart_to_sphe_inv

  • extrapolate_data()

  • find_rotation()

  • s_inv

Calls:

  • dgemm()

  • dgesvd()

get_pseudo_inverse_complex:

File : utils/linear_algebra.irp.f

subroutine get_pseudo_inverse_complex(A,LDA,m,n,C,LDC,cutoff)

Find C = A^-1

Called by:

  • s_inv_complex

Calls:

  • zgesvd()

get_total_available_memory:

File : utils/memory.irp.f

integer function get_total_available_memory() result(res)

Returns the total available memory on the current machine

give_degen:

File : utils/util.irp.f

subroutine give_degen(A, n, shift, list_degen, n_degen_list)

returns n_degen_list :: the number of degenerated SET of elements (i.e. with |A(i)-A(i+1)| below shift)

for each of these sets, list_degen(1,i) = first degenerate element of the set i,

list_degen(2,i) = last degenerate element of the set i.

give_degen_full_list:

File : utils/block_diag_degen.irp.f

subroutine give_degen_full_list(A, n, thr, list_degen, n_degen_list)

you enter with an array A(n) and spits out all the elements degenerated up to thr

the elements of A(n) DON’T HAVE TO BE SORTED IN THE ENTRANCE: TOTALLY GENERAL

list_degen(i,0) = number of degenerate entries

list_degen(i,1) = index of the first degenerate entry

list_degen(i,2:list_degen(i,0)) = list of all other dengenerate entries

if list_degen(i,0) == 1 it means that there is no degeneracy for that element

Called by:

  • diag_mat_per_fock_degen()

give_degen_full_listcore:

File : utils/block_diag_degen_core.irp.f

subroutine give_degen_full_listcore(A, n, listcore, ncore, thr, list_degen, n_degen_list)

you enter with an array A(n) and spits out all the elements degenerated up to thr

the elements of A(n) DON’T HAVE TO BE SORTED IN THE ENTRANCE: TOTALLY GENERAL

list_degen(i,0) = number of degenerate entries

list_degen(i,1) = index of the first degenerate entry

list_degen(i,2:list_degen(i,0)) = list of all other dengenerate entries

if list_degen(i,0) == 1 it means that there is no degeneracy for that element

if list_degen(i,0) >= 1000 it means that it is core orbitals

Called by:

  • diag_mat_per_fock_degen_core()

give_explicit_poly_and_gaussian:

File : utils/integration.irp.f

subroutine give_explicit_poly_and_gaussian(P_new,P_center,p,fact_k,iorder,alpha,beta,a,b,A_center,B_center,dim)
Transforms the product of

(x-x_A)^a(1) (x-x_B)^b(1) (x-x_A)^a(2) (y-y_B)^b(2) (z-z_A)^a(3) (z-z_B)^b(3) exp(-(r-A)^2 alpha) exp(-(r-B)^2 beta)

into
fact_k * [ sum (l_x = 0,i_order(1)) P_new(l_x,1) * (x-P_center(1))^l_x ] exp (- p (x-P_center(1))^2 )
  • [ sum (l_y = 0,i_order(2)) P_new(l_y,2) * (y-P_center(2))^l_y ] exp (- p (y-P_center(2))^2 )

  • [ sum (l_z = 0,i_order(3)) P_new(l_z,3) * (z-P_center(3))^l_z ] exp (- p (z-P_center(3))^2 )

WARNING ::: IF fact_k is too smal then: returns a “s” function centered in zero with an inifinite exponent and a zero polynom coef

Called by:

  • ao_two_e_integral()

  • ao_two_e_integral_erf()

  • ao_two_e_integral_schwartz_accel()

  • ao_two_e_integral_schwartz_accel_erf()

  • give_explicit_poly_and_gaussian_double()

  • overlap_gaussian_xyz()

Calls:

  • gaussian_product()

  • multiply_poly()

  • recentered_poly2()

give_explicit_poly_and_gaussian_double:

File : utils/integration.irp.f

subroutine give_explicit_poly_and_gaussian_double(P_new,P_center,p,fact_k,iorder,alpha,beta,gama,a,b,A_center,B_center,Nucl_center,dim)
Transforms the product of

(x-x_A)^a(1) (x-x_B)^b(1) (x-x_A)^a(2) (y-y_B)^b(2) (z-z_A)^a(3) (z-z_B)^b(3) exp(-(r-A)^2 alpha) exp(-(r-B)^2 beta) exp(-(r-Nucl_center)^2 gama

into
fact_k * [ sum (l_x = 0,i_order(1)) P_new(l_x,1) * (x-P_center(1))^l_x ] exp (- p (x-P_center(1))^2 )
  • [ sum (l_y = 0,i_order(2)) P_new(l_y,2) * (y-P_center(2))^l_y ] exp (- p (y-P_center(2))^2 )

  • [ sum (l_z = 0,i_order(3)) P_new(l_z,3) * (z-P_center(3))^l_z ] exp (- p (z-P_center(3))^2 )

Calls:

  • gaussian_product()

  • give_explicit_poly_and_gaussian()

give_explicit_poly_and_gaussian_v:

File : utils/integration.irp.f

subroutine give_explicit_poly_and_gaussian_v(P_new, ldp, P_center, p, fact_k, iorder, alpha, beta, a, b, A_center, LD_A, B_center, n_points)
Transforms the product of

(x-x_A)^a(1) (x-x_B)^b(1) (y-y_A)^a(2) (y-y_B)^b(2) (z-z_A)^a(3) (z-z_B)^b(3) exp(-(r-A)^2 alpha) exp(-(r-B)^2 beta)

into
fact_k * [ sum (l_x = 0,i_order(1)) P_new(l_x,1) * (x-P_center(1))^l_x ] exp (- p (x-P_center(1))^2 )
  • [ sum (l_y = 0,i_order(2)) P_new(l_y,2) * (y-P_center(2))^l_y ] exp (- p (y-P_center(2))^2 )

  • [ sum (l_z = 0,i_order(3)) P_new(l_z,3) * (z-P_center(3))^l_z ] exp (- p (z-P_center(3))^2 )

WARNING :: : IF fact_k is too smal then: returns a “s” function centered in zero with an inifinite exponent and a zero polynom coef

Called by:

  • overlap_gaussian_xyz_v()

Calls:

  • gaussian_product_v()

  • multiply_poly_v()

  • recentered_poly2_v()

  • recentered_poly2_v0()

give_explicit_poly_and_gaussian_x:

File : utils/integration.irp.f

subroutine give_explicit_poly_and_gaussian_x(P_new,P_center,p,fact_k,iorder,alpha,beta,a,b,A_center,B_center,dim)
Transform the product of

(x-x_A)^a(1) (x-x_B)^b(1) (x-x_A)^a(2) (y-y_B)^b(2) (z-z_A)^a(3) (z-z_B)^b(3) exp(-(r-A)^2 alpha) exp(-(r-B)^2 beta)

into

fact_k (x-x_P)^iorder(1) (y-y_P)^iorder(2) (z-z_P)^iorder(3) exp(-p(r-P)^2)

Called by:

  • overlap_gaussian_x()

Calls:

  • multiply_poly()

  • recentered_poly2()

give_pol_in_r:

File : utils/prim_in_r.irp.f

double precision function give_pol_in_r(r,pol,center, alpha,iorder, max_dim)
hermite:

File : utils/integration.irp.f

double precision function hermite(n,x)

Hermite polynomial

i2set_order:

File : utils/sort.irp.f_template_90

subroutine i2set_order(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A.

i2set_order_big:

File : utils/sort.irp.f_template_90

subroutine i2set_order_big(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A. This is a version for very large arrays where the indices need to be in integer*8 format

i8set_order:

File : utils/sort.irp.f_template_90

subroutine i8set_order(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A.

i8set_order_big:

File : utils/sort.irp.f_template_90

subroutine i8set_order_big(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A. This is a version for very large arrays where the indices need to be in integer*8 format

is_same_spin:

File : utils/util.irp.f

logical function is_same_spin(sigma_1, sigma_2)

true if sgn(sigma_1) = sgn(sigma_2)

iset_order:

File : utils/sort.irp.f_template_90

subroutine iset_order(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A.

Called by:

  • psi_bilinear_matrix_transp_values

  • psi_bilinear_matrix_values

  • restore_symmetry()

iset_order_big:

File : utils/sort.irp.f_template_90

subroutine iset_order_big(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A. This is a version for very large arrays where the indices need to be in integer*8 format

kronecker_delta:

File : utils/util.irp.f

function Kronecker_delta(i, j) result(delta)

Kronecker Delta

lapack_diag:

File : utils/linear_algebra.irp.f

subroutine lapack_diag(eigvalues,eigvectors,H,nmax,n)

Diagonalize matrix H

H is untouched between input and ouptut

eigevalues(i) = ith lowest eigenvalue of the H matrix

eigvectors(i,j) = <i|psi_j> where i is the basis function and psi_j is the j th eigenvector

Called by:

  • ci_electronic_energy

  • davidson_general()

  • davidson_general_ext_rout()

  • davidson_general_ext_rout_diag_dressed()

  • difference_dm_eigvect

  • mo_as_eigvectors_of_mo_matrix()

  • multi_s_x_dipole_moment_eigenvec

  • psi_coef_cas_diagonalized

  • sxeigenvec

Calls:

  • dsyev()

lapack_diag_complex:

File : utils/linear_algebra.irp.f

subroutine lapack_diag_complex(eigvalues,eigvectors,H,nmax,n)

Diagonalize matrix H (complex)

H is untouched between input and ouptut

eigevalues(i) = ith lowest eigenvalue of the H matrix

eigvectors(i,j) = <i|psi_j> where i is the basis function and psi_j is the j th eigenvector

Calls:

  • zheev()

lapack_diagd:

File : utils/linear_algebra.irp.f

subroutine lapack_diagd(eigvalues,eigvectors,H,nmax,n)

Diagonalize matrix H

H is untouched between input and ouptut

eigevalues(i) = ith lowest eigenvalue of the H matrix

eigvectors(i,j) = <i|psi_j> where i is the basis function and psi_j is the j th eigenvector

Called by:

  • eigsvd()

  • exp_matrix()

  • inertia_tensor_eigenvectors

Calls:

  • dsyevd()

lapack_diagd_complex:

File : utils/linear_algebra.irp.f

subroutine lapack_diagd_complex(eigvalues,eigvectors,H,nmax,n)

Diagonalize matrix H(complex)

H is untouched between input and ouptut

eigevalues(i) = ith lowest eigenvalue of the H matrix

eigvectors(i,j) = <i|psi_j> where i is the basis function and psi_j is the j th eigenvector

Calls:

  • zheevd()

lapack_diagd_diag_complex:

File : utils/linear_algebra.irp.f

subroutine lapack_diagd_diag_complex(eigvalues,eigvectors,H,nmax,n)

Diagonalize matrix H(complex)

H is untouched between input and ouptut

eigevalues(i) = ith lowest eigenvalue of the H matrix

eigvectors(i,j) = <i|psi_j> where i is the basis function and psi_j is the j th eigenvector

Calls:

  • zheev()

  • zheevd()

lapack_diagd_diag_in_place_complex:

File : utils/linear_algebra.irp.f

subroutine lapack_diagd_diag_in_place_complex(eigvalues,eigvectors,nmax,n)

Diagonalize matrix H(complex)

H is untouched between input and ouptut

eigevalues(i) = ith lowest eigenvalue of the H matrix

eigvectors(i,j) = <i|psi_j> where i is the basis function and psi_j is the j th eigenvector

Calls:

  • zheev()

  • zheevd()

logfact:

File : utils/util.irp.f

double precision function logfact(n)

n!

lowercase:

File : utils/util.irp.f

subroutine lowercase(txt,n)

Transform to lower case

Called by:

  • end_parallel_job()

  • new_parallel_job()

map_load_from_disk:

File : utils/map_functions.irp.f

subroutine map_load_from_disk(filename,map)

Called by:

  • ao_two_e_integrals_erf_in_map

  • ao_two_e_integrals_in_map

  • mo_two_e_integrals_erf_in_map

  • mo_two_e_integrals_in_map

Calls:

  • c_f_pointer()

  • mmap()

map_save_to_disk:

File : utils/map_functions.irp.f

subroutine map_save_to_disk(filename,map)

Called by:

  • ao_two_e_integrals_erf_in_map

  • ao_two_e_integrals_in_map

  • mo_two_e_integrals_erf_in_map

  • mo_two_e_integrals_in_map

  • save_erf_two_e_integrals_ao()

  • save_erf_two_e_integrals_mo()

  • save_erf_two_e_ints_ao_into_ints_ao()

  • save_erf_two_e_ints_mo_into_ints_mo()

Calls:

  • c_f_pointer()

  • map_sort()

  • mmap()

  • msync()

matrix_vector_product_complex:

File : utils/linear_algebra.irp.f

subroutine matrix_vector_product_complex(u0,u1,matrix,sze,lda)

performs u1 =! performs u1 +( u0 * matrix)

Calls:

  • zhemv()

memory_of_double:

File : utils/memory.irp.f

double precision function memory_of_double(n)

Computes the memory required for n double precision elements in gigabytes.

memory_of_double8:

File : utils/memory.irp.f

double precision function memory_of_double8(n)

Computes the memory required for n double precision elements in gigabytes.

memory_of_int:

File : utils/memory.irp.f

double precision function memory_of_int(n)

Computes the memory required for n double precision elements in gigabytes.

memory_of_int8:

File : utils/memory.irp.f

double precision function memory_of_int8(n)

Computes the memory required for n double precision elements in gigabytes.

multiply_cpoly:

File : utils/cgtos_utils.irp.f

subroutine multiply_cpoly(b, nb, c, nc, d, nd)

Multiply two complex polynomials D(t) =! D(t) +( B(t) * C(t))

Called by:

  • general_primitive_integral_cgtos()

  • give_cpolynomial_mult_center_one_e()

  • give_explicit_cpoly_and_cgaussian()

  • give_explicit_cpoly_and_cgaussian_x()

  • i_x1_pol_mult_a1_cgtos()

  • i_x1_pol_mult_a2_cgtos()

  • i_x1_pol_mult_one_e_cgtos()

  • i_x1_pol_mult_recurs_cgtos()

  • i_x2_pol_mult_cgtos()

  • i_x2_pol_mult_one_e_cgtos()

multiply_poly:

File : utils/integration.irp.f

subroutine multiply_poly(b,nb,c,nc,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • general_primitive_integral_erf()

  • give_explicit_poly_and_gaussian()

  • give_explicit_poly_and_gaussian_x()

  • give_polynomial_mult_center_one_e()

  • give_polynomial_mult_center_one_e_erf()

  • give_polynomial_mult_center_one_e_erf_opt()

Calls:

  • multiply_poly_b0()

  • multiply_poly_b1()

  • multiply_poly_b2()

  • multiply_poly_c0()

  • multiply_poly_c1()

  • multiply_poly_c2()

multiply_poly_b0:

File : utils/integration.irp.f

subroutine multiply_poly_b0(b,c,nc,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • multiply_poly()

multiply_poly_b1:

File : utils/integration.irp.f

subroutine multiply_poly_b1(b,c,nc,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • multiply_poly()

multiply_poly_b2:

File : utils/integration.irp.f

subroutine multiply_poly_b2(b,c,nc,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • multiply_poly()

multiply_poly_c0:

File : utils/integration.irp.f

subroutine multiply_poly_c0(b,nb,c,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • multiply_poly()

multiply_poly_c1:

File : utils/integration.irp.f

subroutine multiply_poly_c1(b,nb,c,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • multiply_poly()

multiply_poly_c2:

File : utils/integration.irp.f

subroutine multiply_poly_c2(b,nb,c,d,nd)

Multiply two polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • i_x1_pol_mult_one_e()

  • i_x2_pol_mult_one_e()

  • multiply_poly()

multiply_poly_v:

File : utils/integration.irp.f

subroutine multiply_poly_v(b,nb,c,nc,d,nd,n_points)

Multiply pairs of polynomials D(t) =! D(t) +( B(t)*C(t))

Called by:

  • give_explicit_poly_and_gaussian_v()

normalize:

File : utils/util.irp.f

subroutine normalize(u,sze)

Normalizes vector u

Called by:

  • copy_h_apply_buffer_to_wf()

  • dav_double_dressed()

  • davidson_diag_csf_hjj()

  • davidson_diag_hjj()

  • davidson_diag_hjj_sjj()

  • davidson_diag_nonsym_hjj()

  • davidson_general()

  • davidson_general_diag_dressed_ext_rout_nonsym_b1space()

  • davidson_general_ext_rout()

  • davidson_general_ext_rout_diag_dressed()

  • davidson_general_ext_rout_dressed()

  • davidson_general_ext_rout_nonsym_b1space()

  • save_wavefunction_general()

Calls:

  • dscal()

nullify_small_elements:

File : utils/linear_algebra.irp.f

subroutine nullify_small_elements(m,n,A,LDA,thresh)

Called by:

  • ci_electronic_energy

  • dav_double_dressed()

  • davidson_diag_csf_hjj()

  • davidson_diag_hjj()

  • davidson_diag_hjj_sjj()

  • davidson_diag_nonsym_hjj()

  • davidson_general_ext_rout_dressed()

  • hcore_guess()

  • mo_one_e_integrals

  • save_natural_mos()

  • save_natural_mos_canon_label()

  • save_natural_mos_no_ov_rot()

  • save_wavefunction_truncated()

ortho_canonical:

File : utils/linear_algebra.irp.f

subroutine ortho_canonical(overlap,LDA,N,C,LDC,m,cutoff)

Compute C_new=C_old.U.s^-1/2 canonical orthogonalization.

overlap : overlap matrix

LDA : leftmost dimension of overlap array

N : Overlap matrix is NxN (array is (LDA,N) )

CCoefficients of the vectors to orthogonalize. On exit,

orthogonal vectors

LDC : leftmost dimension of C

m : Coefficients matrix is MxN, ( array is (LDC,N) )

Called by:

  • ao_ortho_canonical_coef

Calls:

  • dgemm()

  • svd()

ortho_canonical_complex:

File : utils/linear_algebra.irp.f

subroutine ortho_canonical_complex(overlap,LDA,N,C,LDC,m,cutoff)

Compute C_new=C_old.U.s^-1/2 canonical orthogonalization.

overlap : overlap matrix

LDA : leftmost dimension of overlap array

N : Overlap matrix is NxN (array is (LDA,N) )

CCoefficients of the vectors to orthogonalize. On exit,

orthogonal vectors

LDC : leftmost dimension of C

m : Coefficients matrix is MxN, ( array is (LDC,N) )

Calls:

  • svd_complex()

  • zgemm()

ortho_lowdin:

File : utils/linear_algebra.irp.f

subroutine ortho_lowdin(overlap,LDA,N,C,LDC,m,cutoff)

Compute C_new=C_old.S^-1/2 orthogonalization.

overlap : overlap matrix

LDA : leftmost dimension of overlap array

N : Overlap matrix is NxN (array is (LDA,N) )

CCoefficients of the vectors to orthogonalize. On exit,

orthogonal vectors

LDC : leftmost dimension of C

M : Coefficients matrix is MxN, ( array is (LDC,N) )

Called by:

  • ao_ortho_lowdin_coef

  • orthonormalize_mos()

Calls:

  • dgemm()

  • svd()

ortho_lowdin_complex:

File : utils/linear_algebra.irp.f

subroutine ortho_lowdin_complex(overlap,LDA,N,C,LDC,m,cutoff)

Compute C_new=C_old.S^-1/2 orthogonalization.

overlap : overlap matrix

LDA : leftmost dimension of overlap array

N : Overlap matrix is NxN (array is (LDA,N) )

CCoefficients of the vectors to orthogonalize. On exit,

orthogonal vectors

LDC : leftmost dimension of C

M : Coefficients matrix is MxN, ( array is (LDC,N) )

Calls:

  • svd_complex()

  • zgemm()

ortho_qr:

File : utils/linear_algebra.irp.f

subroutine ortho_qr(A,LDA,m,n)

Orthogonalization using Q.R factorization

A : matrix to orthogonalize

LDA : leftmost dimension of A

m : Number of rows of A

n : Number of columns of A

Called by:

  • davidson_diag_nonsym_hjj()

  • davidson_general()

  • davidson_general_diag_dressed_ext_rout_nonsym_b1space()

  • davidson_general_ext_rout()

  • davidson_general_ext_rout_diag_dressed()

  • davidson_general_ext_rout_dressed()

  • davidson_general_ext_rout_nonsym_b1space()

  • ortho_svd()

Calls:

  • dgeqrf()

  • dorgqr()

ortho_qr_complex:

File : utils/linear_algebra.irp.f

subroutine ortho_qr_complex(A,LDA,m,n)

Orthogonalization using Q.R factorization

A : matrix to orthogonalize

LDA : leftmost dimension of A

n : Number of rows of A

m : Number of columns of A

Calls:

  • zgeqrf()

  • zungqr()

ortho_qr_unblocked:

File : utils/linear_algebra.irp.f

subroutine ortho_qr_unblocked(A,LDA,m,n)

Orthogonalization using Q.R factorization

A : matrix to orthogonalize

LDA : leftmost dimension of A

n : Number of rows of A

m : Number of columns of A

Calls:

  • dgeqr2()

  • dorg2r()

ortho_qr_unblocked_complex:

File : utils/linear_algebra.irp.f

subroutine ortho_qr_unblocked_complex(A,LDA,m,n)

Orthogonalization using Q.R factorization

A : matrix to orthogonalize

LDA : leftmost dimension of A

n : Number of rows of A

m : Number of columns of A

ortho_svd:

File : utils/linear_algebra.irp.f

subroutine ortho_svd(A,LDA,m,n)

Orthogonalization via fast SVD

A : matrix to orthogonalize

LDA : leftmost dimension of A

m : Number of rows of A

n : Number of columns of A

Called by:

  • randomized_svd()

Calls:

  • ortho_qr()

  • svd()

overlap_cgaussian_x:

File : utils/cgtos_one_e.irp.f

complex*16 function overlap_cgaussian_x(Ae_center, Be_center, alpha, beta, power_A, power_B, Ap_center, Bp_center, dim)

int_{-infty}^{+infty} (x - Ap_x)^ax (x - Bp_x)^bx exp(-alpha (x - Ae_x)^2) exp(-beta (x - Be_X)^2) dx

with complex arguments

Calls:

  • give_explicit_cpoly_and_cgaussian_x()

overlap_gaussian_x:

File : utils/one_e_integration.irp.f

double precision function overlap_gaussian_x(A_center,B_center,alpha,beta,power_A,power_B,dim)
\[\sum_{-infty}^{+infty} (x-A_x)^ax (x-B_x)^bx exp(-alpha(x-A_x)^2) exp(-beta(x-B_X)^2) dx\]

Calls:

  • give_explicit_poly_and_gaussian_x()

overlap_gaussian_xyz:

File : utils/one_e_integration.irp.f

subroutine overlap_gaussian_xyz(A_center, B_center, alpha, beta, power_A, power_B, overlap_x, overlap_y, overlap_z, overlap, dim)
\[\begin{split}S_x = \int (x-A_x)^{a_x} exp(-\alpha(x-A_x)^2) (x-B_x)^{b_x} exp(-beta(x-B_x)^2) dx \\ S = S_x S_y S_z\end{split}\]

Called by:

  • ao_coef_normalized

  • ao_deriv2_x

  • ao_deriv_1_x

  • ao_dipole_x

  • ao_overlap

  • ao_spread_x

  • prim_normalization_factor

  • shell_normalization_factor

Calls:

  • gaussian_product_x()

  • give_explicit_poly_and_gaussian()

overlap_gaussian_xyz_v:

File : utils/one_e_integration.irp.f

subroutine overlap_gaussian_xyz_v(A_center, B_center, alpha, beta, power_A, power_B, overlap, n_points)
\[\begin{split}S_x = \int (x-A_x)^{a_x} exp(-\alpha(x-A_x)^2) (x-B_x)^{b_x} exp(-beta(x-B_x)^2) dx \\ S = S_x S_y S_z\end{split}\]

Calls:

  • give_explicit_poly_and_gaussian_v()

overlap_x_abs:

File : utils/one_e_integration.irp.f

subroutine overlap_x_abs(A_center, B_center, alpha, beta, power_A, power_B, overlap_x, lower_exp_val, dx, nx)

Called by:

  • ao_overlap_abs

pivoted_cholesky:

File : utils/linear_algebra.irp.f

subroutine pivoted_cholesky( A, rank, tol, ndim, U)

Called by:

  • roothaan_hall_scf()

Calls:

  • dpstrf()

pol_modif_center:

File : utils/integration.irp.f

subroutine pol_modif_center(A_center, B_center, iorder, A_pol, B_pol)
Transform the pol centerd on A:

[ sum_i ax_i (x-x_A)^i ] [ sum_j ay_j (y-y_A)^j ] [ sum_k az_k (z-z_A)^k ]

to a pol centered on B

[ sum_i bx_i (x-x_B)^i ] [ sum_j by_j (y-y_B)^j ] [ sum_k bz_k (z-z_B)^k ]

Calls:

  • pol_modif_center_x()

pol_modif_center_x:

File : utils/integration.irp.f

subroutine pol_modif_center_x(A_center, B_center, iorder, A_pol, B_pol)
Transform the pol centerd on A:

[ sum_i ax_i (x-x_A)^i ]

to a pol centered on B

[ sum_i bx_i (x-x_B)^i ]

bx_i = sum_{j=i}^{iorder} ax_j (x_B - x_A)^(j-i) j! / [ i! (j-i)! ]

= sum_{j=i}^{iorder} ax_j (x_B - x_A)^(j-i) binom_func(j,i)

Called by:

  • pol_modif_center()

primitive_value_explicit:

File : utils/prim_in_r.irp.f

double precision function primitive_value_explicit(power_prim,center_prim,alpha,r)

Evaluates at “r” a primitive of type : (x - center_prim(1))**power_prim(1) (y - center_prim(2))**power_prim(2) * (z - center_prim(3))**power_prim(3)

exp(-alpha * [(x - center_prim(1))**2 + (y - center_prim(2))**2 + (z - center_prim(3))**2] )

print_memory_usage:

File : utils/memory.irp.f

subroutine print_memory_usage()

Prints the memory usage in the output

Called by:

  • check_mem()

  • cholesky_ao_num

  • write_time()

Calls:

  • resident_memory()

  • total_memory()

randomized_svd:

File : utils/linear_algebra.irp.f

subroutine randomized_svd(A,LDA,U,LDU,D,Vt,LDVt,m,n,q,r)

Randomized SVD: rank r, q power iterations

  1. Sample column space of A with P: Z = A.P where P is random with r+p columns.

  2. Power iterations : Z <- X . (Xt.Z)

  3. Z = Q.R

  4. Compute SVD on projected Qt.X = U’ . S. Vt

  5. U = Q U’

Calls:

  • dgemm()

  • ortho_svd()

  • random_number()

  • svd()

recentered_cpoly2:

File : utils/cgtos_utils.irp.f

subroutine recentered_cpoly2(P_A, x_A, x_P, a, P_B, x_B, x_Q, b)

write two complex polynomials (x-x_A)^a (x-x_B)^b as P_A(x-x_P) and P_B(x-x_Q)

Needs:

  • binom

Called by:

  • give_explicit_cpoly_and_cgaussian()

  • give_explicit_cpoly_and_cgaussian_x()

recentered_poly2:

File : utils/integration.irp.f

subroutine recentered_poly2(P_new, x_A, x_P, a, Q_new, x_B, x_Q, b)

Recenter two polynomials:

(x - x_A)^a -> sum_{i=0}^{a} binom{a}{i} (x_A - x_P)^{a-i} (x - x_P)^i (x - x_B)^b -> sum_{i=0}^{b} binom{b}{i} (x_B - x_Q)^{b-i} (x - x_Q)^i

where:

binom{a}{i} = binom{a}{a-i} = a! / [i! (a-i)!]

Needs:

  • binom

Called by:

  • give_explicit_poly_and_gaussian()

  • give_explicit_poly_and_gaussian_x()

recentered_poly2_v:

File : utils/integration.irp.f

subroutine recentered_poly2_v(P_new, lda, x_A, LD_xA, x_P, a, P_new2, ldb, x_B, x_Q, b, n_points)

Recenter two polynomials

Needs:

  • binom

Called by:

  • give_explicit_poly_and_gaussian_v()

recentered_poly2_v0:

File : utils/integration.irp.f

subroutine recentered_poly2_v0(P_new, lda, x_A, LD_xA, x_P, a, n_points)

Recenter two polynomials. Special case for b=(0,0,0)

(x - A)^a (x - B)^0 = (x - P + P - A)^a (x - Q + Q - B)^0

= (x - P + P - A)^a

Needs:

  • binom

Called by:

  • give_explicit_poly_and_gaussian_v()

resident_memory:

File : utils/memory.irp.f

subroutine resident_memory(value)

Returns the current used memory in gigabytes used by the current process.

Called by:

  • ao_two_e_integral_alpha_chol

  • check_mem()

  • cholesky_ao_num

  • dav_double_dressed()

  • davidson_diag_csf_hjj()

  • davidson_diag_hjj()

  • davidson_diag_hjj_sjj()

  • davidson_diag_nonsym_hjj()

  • davidson_general()

  • davidson_general_diag_dressed_ext_rout_nonsym_b1space()

  • davidson_general_ext_rout()

  • davidson_general_ext_rout_diag_dressed()

  • davidson_general_ext_rout_dressed()

  • davidson_general_ext_rout_nonsym_b1space()

  • print_memory_usage()

  • run_slave_main()

  • zmq_pt2()

Calls:

  • lock_io()

  • unlock_io()

  • usleep()

restore_symmetry:

File : utils/linear_algebra.irp.f

subroutine restore_symmetry(m,n,A,LDA,thresh)

Tries to find the matrix elements that are the same, and sets them to the average value. If restore_symm is False, only nullify small elements

Called by:

  • ao_to_mo()

  • create_guess()

  • huckel_guess()

  • roothaan_hall_scf()

  • svd_symm()

Calls:

  • dset_order()

  • dsort()

  • iset_order()

rint:

File : utils/integration.irp.f

double precision function rint(n,rho)
\[\int_0^1 dx \exp(-p x^2) x^n\]
rint1:

File : utils/integration.irp.f

double precision function rint1(n,rho)

Standard version of rint

Needs:

  • fact_inv

  • inv_int

rint_large_n:

File : utils/integration.irp.f

double precision function rint_large_n(n,rho)

Version of rint for large values of n

rint_sum:

File : utils/integration.irp.f

double precision function rint_sum(n_pt_out,rho,d1)

Needed for the calculation of two-electron integrals.

set_multiple_levels_omp:

File : utils/set_multiple_levels_omp.irp.f

subroutine set_multiple_levels_omp(activate)

If true, activate OpenMP nested parallelism. If false, deactivate.

Called by:

  • add_integrals_to_map_cholesky()

  • cholesky_ao_num

  • cholesky_mo

  • h_s2_u_0_nstates_zmq()

  • h_u_0_nstates_zmq()

  • mo_integrals_cache

  • run_slave_cipsi()

  • run_slave_main()

  • zmq_pt2()

Calls:

  • omp_set_max_active_levels()

  • omp_set_nested()

set_order:

File : utils/sort.irp.f_template_90

subroutine set_order(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A.

set_order_big:

File : utils/sort.irp.f_template_90

subroutine set_order_big(x,iorder,isize)

array A has already been sorted, and iorder has contains the new order of elements of A. This subroutine changes the order of x to match the new order of A. This is a version for very large arrays where the indices need to be in integer*8 format

shank:

File : utils/shank.irp.f

subroutine shank(array,n,nmax,shank1)

Called by:

  • shank_general()

shank_function:

File : utils/shank.irp.f

double precision function shank_function(array,i,n)
shank_general:

File : utils/shank.irp.f

double precision function shank_general(array,n,nmax)

Calls:

  • shank()

sorted_dnumber:

File : utils/sort.irp.f_template_32

subroutine sorted_dnumber(x,isize,n)

Returns the number of sorted elements

sorted_i2number:

File : utils/sort.irp.f_template_32

subroutine sorted_i2number(x,isize,n)

Returns the number of sorted elements

sorted_i8number:

File : utils/sort.irp.f_template_32

subroutine sorted_i8number(x,isize,n)

Returns the number of sorted elements

sorted_inumber:

File : utils/sort.irp.f_template_32

subroutine sorted_inumber(x,isize,n)

Returns the number of sorted elements

sorted_number:

File : utils/sort.irp.f_template_32

subroutine sorted_number(x,isize,n)

Returns the number of sorted elements

sub_a_at:

File : utils/util.irp.f

subroutine sub_A_At(A, N)
sum_a_at:

File : utils/util.irp.f

subroutine sum_A_At(A, N)
svd:

File : utils/linear_algebra.irp.f

subroutine svd(A,LDA,U,LDU,D,Vt,LDVt,m,n)

Compute A = U.D.Vt

LDx : leftmost dimension of x

Dimension of A is m x n

Called by:

  • eigsvd()

  • mo_as_svd_vectors_of_mo_matrix()

  • mo_as_svd_vectors_of_mo_matrix_eig()

  • mo_coef_new_as_svd_vectors_of_mo_matrix_eig()

  • natorbsci

  • ortho_canonical()

  • ortho_lowdin()

  • ortho_svd()

  • randomized_svd()

  • s_half

  • s_half_inv

Calls:

  • dgesvd()

svd_complex:

File : utils/linear_algebra.irp.f

subroutine svd_complex(A,LDA,U,LDU,D,Vt,LDVt,m,n)

Compute A = U.D.Vt

LDx : leftmost dimension of x

Dimension of A is m x n A,U,Vt are complex*16 D is double precision

Called by:

  • ortho_canonical_complex()

  • ortho_lowdin_complex()

Calls:

  • zgesvd()

svd_symm:

File : utils/linear_algebra.irp.f

subroutine svd_symm(A,LDA,U,LDU,D,Vt,LDVt,m,n)

Compute A = U.D.Vt

LDx : leftmost dimension of x

Dimension of A is m x n

Calls:

  • dgemm()

  • dgesvd()

  • restore_symmetry()

total_memory:

File : utils/memory.irp.f

subroutine total_memory(value)

Returns the current used memory in gigabytes used by the current process.

Called by:

  • print_memory_usage()

Calls:

  • lock_io()

  • unlock_io()

u_dot_u:

File : utils/util.irp.f

double precision function u_dot_u(u,sze)

Compute <u|u>

u_dot_v:

File : utils/util.irp.f

double precision function u_dot_v(u,v,sze)

Compute <u|v>

v2_over_x:

File : utils/util.irp.f

subroutine v2_over_x(v,x,res)
v_phi:

File : utils/integration.irp.f

double precision function V_phi(n, m)

Computes the angular $phi$ part of the nuclear attraction integral:

$int_{0}^{2 pi} cos(phi)^n sin(phi)^m dphi$.

v_theta:

File : utils/integration.irp.f

double precision function V_theta(n, m)

Computes the angular $theta$ part of the nuclear attraction integral:

$int_{0}^{pi} cos(theta)^n sin(theta)^m dtheta$

wall_time:

File : utils/util.irp.f

subroutine wall_time(t)

The equivalent of cpu_time, but for the wall time.

Called by:

  • act_2_rdm_aa_mo

  • act_2_rdm_ab_mo

  • act_2_rdm_bb_mo

  • act_2_rdm_spin_trace_mo

  • act_2_rdm_trans_spin_trace_mo

  • add_integrals_to_map()

  • add_integrals_to_map_erf()

  • ao_pseudo_integrals_local

  • ao_pseudo_integrals_non_local

  • ao_two_e_integrals_erf_in_map

  • ao_two_e_integrals_in_map

  • apply_mo_rotation()

  • bielecci

  • bielecci_no

  • cholesky_ao_num

  • cholesky_mo_transp

  • cholesky_no_1_idx_transp

  • cholesky_no_2_idx_transp

  • cholesky_no_total_transp

  • cholesky_semi_mo_transp_simple

  • diag_hessian_list_opt()

  • diag_hessian_opt()

  • diagonalization_hessian()

  • first_diag_hessian_list_opt()

  • first_diag_hessian_opt()

  • first_hessian_list_opt()

  • first_hessian_opt()

  • gradient_list_opt()

  • gradient_opt()

  • hessian_list_opt()

  • hessian_opt()

  • mo_two_e_integrals_erf_in_map

  • mo_two_e_integrals_in_map

  • output_wall_time_0

  • pt2_collector()

  • rotation_matrix()

  • rotation_matrix_iterative()

  • run_pt2_slave_large()

  • run_pt2_slave_small()

  • run_slave_main()

  • state_av_act_2_rdm_aa_mo

  • state_av_act_2_rdm_ab_mo

  • state_av_act_2_rdm_bb_mo

  • state_av_act_2_rdm_spin_trace_mo

  • trust_region_expected_e()

  • trust_region_optimal_lambda()

  • trust_region_step()

  • write_time()

Calls:

  • system_clock()

wallis:

File : utils/integration.irp.f

double precision function Wallis(n)

Wallis integral:

$int_{0}^{pi} cos(theta)^n dtheta$.

write_git_log:

File : utils/util.irp.f

subroutine write_git_log(iunit)

Write the last git commit in file iunit.

zboysfun:

File : utils/cpx_boys.irp.f

subroutine zboysfun(n_max, x, vals)

Computes values of the Boys function for n = 0, 1, …, n_max for a complex valued argument

Input: x — argument, complex*16, Re(x) >= 0 Output: vals — values of the Boys function, n = 0, 1, …, n_max

Beylkin & Sharma, J. Chem. Phys. 155, 174117 (2021) https://doi.org/10.1063/5.0062444

Called by:

  • crint_2()

  • crint_2_vec()

Calls:

  • zboysfun00_2()

zboysfun00_1:

File : utils/cpx_erf.irp.f

subroutine zboysfun00_1(rho, F0)

Called by:

  • crint_1()

  • crint_1_vec()

zboysfun00_2:

File : utils/cpx_erf.irp.f

subroutine zboysfun00_2(z, val)

Computes values of the Boys function for n=0 for a complex valued argument

Input: z — argument, complex*16, Real(z) >= 0 Output: val — value of the Boys function n=0

Beylkin & Sharma, J. Chem. Phys. 155, 174117 (2021) https://doi.org/10.1063/5.0062444

Called by:

  • zboysfun()

zboysfun00nrp:

File : utils/cpx_erf.irp.f

subroutine zboysfun00nrp(z, val)

Computes values of the exp(z) F(0,z) (where F(0,z) is the Boys function) for a complex valued argument with Real(z)<=0

Input: z — argument, complex*16, !!! Real(z)<=0 !!! Output: val — value of the function !!! exp(z) F(0,z) !!!, where F(0,z) is the Boys function

Beylkin & Sharma, J. Chem. Phys. 155, 174117 (2021) https://doi.org/10.1063/5.0062444

Called by:

  • zboysfunnrp()

zboysfunnrp:

File : utils/cpx_boys.irp.f

subroutine zboysfunnrp(n_max, x, vals)

Computes values of e^z F(n,z) for n = 0, 1, …, n_max (where F(n,z) are the Boys functions) for a complex valued argument WITH NEGATIVE REAL PART

Input: x — argument, complex *16 Re(x)<=0 Output: vals — values of e^z F(n,z), n = 0, 1, …, n_max

Beylkin & Sharma, J. Chem. Phys. 155, 174117 (2021) https://doi.org/10.1063/5.0062444

Called by:

  • crint_2()

  • crint_2_vec()

Calls:

  • zboysfun00nrp()