dft_utils_one_e

This module contains all the one-body related quantities needed to perform DFT or RS-DFT calculations with the LDA and PBE functionals. Therefore, it contains most of the properties which depends on the one-body density and density matrix.

Some interesting quantities you might take a look at:

  • The LDA and PBE providers for the x/c energies in e_xc.irp.f and sr_exc.irp.f

  • The LDA and PBE providers for the x/c potentials on the AO basis in pot_ao.irp.f and sr_pot_ao.irp.f

  • The \(h_{core}\) energy computed directly with the one-body density matrix in one_e_energy_dft.irp.f

  • LDA and PBE short-range functionals subroutines in exc_sr_lda.irp.f and exc_sr_pbe.irp.f

Providers

ao_effective_one_e_potential

File : dft_utils_one_e/effective_pot.irp.f

double precision, allocatable   :: ao_effective_one_e_potential (ao_num,ao_num,N_states)
double precision, allocatable   :: ao_effective_one_e_potential_without_kin     (ao_num,ao_num,N_states)

ao_effective_one_e_potential(i,j) = \(\rangle i_{AO}| v_{H}^{sr} |j_{AO}\rangle + \rangle i_{AO}| h_{core} |j_{AO}\rangle + \rangle i_{AO}|v_{xc} |j_{AO}\rangle\)

Needs:

  • ao_num

  • effective_one_e_potential

  • mo_coef

  • mo_num

  • n_states

ao_effective_one_e_potential_without_kin

File : dft_utils_one_e/effective_pot.irp.f

double precision, allocatable   :: ao_effective_one_e_potential (ao_num,ao_num,N_states)
double precision, allocatable   :: ao_effective_one_e_potential_without_kin     (ao_num,ao_num,N_states)

ao_effective_one_e_potential(i,j) = \(\rangle i_{AO}| v_{H}^{sr} |j_{AO}\rangle + \rangle i_{AO}| h_{core} |j_{AO}\rangle + \rangle i_{AO}|v_{xc} |j_{AO}\rangle\)

Needs:

  • ao_num

  • effective_one_e_potential

  • mo_coef

  • mo_num

  • n_states

effective_one_e_potential

File : dft_utils_one_e/effective_pot.irp.f

double precision, allocatable   :: effective_one_e_potential    (mo_num,mo_num,N_states)
double precision, allocatable   :: effective_one_e_potential_without_kin        (mo_num,mo_num,N_states)

Effective_one_e_potential(i,j) = \(\rangle i_{MO}| v_{H}^{sr} |j_{MO}\rangle + \rangle i_{MO}| h_{core} |j_{MO}\rangle + \rangle i_{MO}|v_{xc} |j_{MO}\rangle\)

on the MO basis Taking the expectation value does not provide any energy, but effective_one_e_potential(i,j) is the potential coupling DFT and WFT part to be used in any WFT calculation.

Needs:

  • mo_integrals_n_e

  • mo_kinetic_integrals

  • mo_num

  • n_states

  • potential_c_alpha_mo

  • potential_x_alpha_mo

  • short_range_hartree_operator

Needed by:

  • ao_effective_one_e_potential

effective_one_e_potential_without_kin

File : dft_utils_one_e/effective_pot.irp.f

double precision, allocatable   :: effective_one_e_potential    (mo_num,mo_num,N_states)
double precision, allocatable   :: effective_one_e_potential_without_kin        (mo_num,mo_num,N_states)

Effective_one_e_potential(i,j) = \(\rangle i_{MO}| v_{H}^{sr} |j_{MO}\rangle + \rangle i_{MO}| h_{core} |j_{MO}\rangle + \rangle i_{MO}|v_{xc} |j_{MO}\rangle\)

on the MO basis Taking the expectation value does not provide any energy, but effective_one_e_potential(i,j) is the potential coupling DFT and WFT part to be used in any WFT calculation.

Needs:

  • mo_integrals_n_e

  • mo_kinetic_integrals

  • mo_num

  • n_states

  • potential_c_alpha_mo

  • potential_x_alpha_mo

  • short_range_hartree_operator

Needed by:

  • ao_effective_one_e_potential

energy_sr_c_lda

File : dft_utils_one_e/sr_exc.irp.f

double precision, allocatable   :: energy_sr_x_lda      (N_states)
double precision, allocatable   :: energy_sr_c_lda      (N_states)

exchange/correlation energy with the short range lda functional

Needs:

  • final_grid_points

  • mu_erf_dft

  • n_points_final_grid

  • n_states

  • one_e_dm_alpha_at_r

energy_sr_c_pbe

File : dft_utils_one_e/sr_exc.irp.f

double precision, allocatable   :: energy_sr_x_pbe      (N_states)
double precision, allocatable   :: energy_sr_c_pbe      (N_states)

exchange/correlation energy with the short range pbe functional

Needs:

  • final_grid_points

  • mu_erf_dft

  • n_points_final_grid

  • n_states

  • one_e_dm_and_grad_alpha_in_r

energy_sr_x_lda

File : dft_utils_one_e/sr_exc.irp.f

double precision, allocatable   :: energy_sr_x_lda      (N_states)
double precision, allocatable   :: energy_sr_c_lda      (N_states)

exchange/correlation energy with the short range lda functional

Needs:

  • final_grid_points

  • mu_erf_dft

  • n_points_final_grid

  • n_states

  • one_e_dm_alpha_at_r

energy_sr_x_pbe

File : dft_utils_one_e/sr_exc.irp.f

double precision, allocatable   :: energy_sr_x_pbe      (N_states)
double precision, allocatable   :: energy_sr_c_pbe      (N_states)

exchange/correlation energy with the short range pbe functional

Needs:

  • final_grid_points

  • mu_erf_dft

  • n_points_final_grid

  • n_states

  • one_e_dm_and_grad_alpha_in_r

gga_sr_type_functionals:

File : dft_utils_one_e/utils.irp.f

subroutine GGA_sr_type_functionals(r,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,grad_rho_a_b, &
                        ex,vx_rho_a,vx_rho_b,vx_grad_rho_a_2,vx_grad_rho_b_2,vx_grad_rho_a_b, &
                        ec,vc_rho_a,vc_rho_b,vc_grad_rho_a_2,vc_grad_rho_b_2,vc_grad_rho_a_b )

routine that helps in building the x/c potentials on the AO basis for a GGA functional with a short-range interaction

Needs:

  • mu_erf_dft

  • n_states

Called by:

  • aos_sr_vc_alpha_pbe_w

  • aos_sr_vxc_alpha_pbe_w

  • energy_sr_x_pbe

  • energy_x_sr_pbe

Calls:

  • ec_pbe_sr()

  • ex_pbe_sr()

  • grad_rho_ab_to_grad_rho_oc()

  • rho_ab_to_rho_oc()

  • v_grad_rho_oc_to_v_grad_rho_ab()

  • v_rho_oc_to_v_rho_ab()

gga_type_functionals:

File : dft_utils_one_e/utils.irp.f

subroutine GGA_type_functionals(r,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,grad_rho_a_b, &
                        ex,vx_rho_a,vx_rho_b,vx_grad_rho_a_2,vx_grad_rho_b_2,vx_grad_rho_a_b, &
                        ec,vc_rho_a,vc_rho_b,vc_grad_rho_a_2,vc_grad_rho_b_2,vc_grad_rho_a_b )

routine that helps in building the x/c potentials on the AO basis for a GGA functional

Needs:

  • n_states

Called by:

  • aos_vc_alpha_pbe_w

  • aos_vxc_alpha_pbe_w

  • energy_c_pbe

  • energy_x_pbe

Calls:

  • ec_pbe_sr()

  • ex_pbe_sr()

  • grad_rho_ab_to_grad_rho_oc()

  • rho_ab_to_rho_oc()

  • v_grad_rho_oc_to_v_grad_rho_ab()

  • v_rho_oc_to_v_rho_ab()

mu_erf_dft

File : dft_utils_one_e/mu_erf_dft.irp.f

double precision        :: mu_erf_dft

range separation parameter used in RS-DFT. It is set to mu_erf in order to be consistent with the two electrons integrals erf

Needs:

  • mu_erf

Needed by:

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vc_alpha_pbe_w

  • aos_sr_vxc_alpha_lda_w

  • aos_sr_vxc_alpha_pbe_w

  • energy_c_sr_lda

  • energy_sr_x_lda

  • energy_sr_x_pbe

  • energy_x_sr_lda

  • energy_x_sr_pbe

psi_dft_energy_h_core

File : dft_utils_one_e/one_e_energy_dft.irp.f

double precision, allocatable   :: psi_dft_energy_kinetic       (N_states)
double precision, allocatable   :: psi_dft_energy_nuclear_elec  (N_states)
double precision, allocatable   :: psi_dft_energy_h_core        (N_states)

kinetic, electron-nuclear and total h_core energy computed with the density matrix one_e_dm_mo_beta_for_dft+one_e_dm_mo_alpha_for_dft

Needs:

  • elec_alpha_num

  • elec_beta_num

  • mo_integrals_n_e

  • mo_kinetic_integrals

  • mo_num

  • n_states

  • one_e_dm_mo_alpha_for_dft

  • one_e_dm_mo_beta_for_dft

psi_dft_energy_kinetic

File : dft_utils_one_e/one_e_energy_dft.irp.f

double precision, allocatable   :: psi_dft_energy_kinetic       (N_states)
double precision, allocatable   :: psi_dft_energy_nuclear_elec  (N_states)
double precision, allocatable   :: psi_dft_energy_h_core        (N_states)

kinetic, electron-nuclear and total h_core energy computed with the density matrix one_e_dm_mo_beta_for_dft+one_e_dm_mo_alpha_for_dft

Needs:

  • elec_alpha_num

  • elec_beta_num

  • mo_integrals_n_e

  • mo_kinetic_integrals

  • mo_num

  • n_states

  • one_e_dm_mo_alpha_for_dft

  • one_e_dm_mo_beta_for_dft

psi_dft_energy_nuclear_elec

File : dft_utils_one_e/one_e_energy_dft.irp.f

double precision, allocatable   :: psi_dft_energy_kinetic       (N_states)
double precision, allocatable   :: psi_dft_energy_nuclear_elec  (N_states)
double precision, allocatable   :: psi_dft_energy_h_core        (N_states)

kinetic, electron-nuclear and total h_core energy computed with the density matrix one_e_dm_mo_beta_for_dft+one_e_dm_mo_alpha_for_dft

Needs:

  • elec_alpha_num

  • elec_beta_num

  • mo_integrals_n_e

  • mo_kinetic_integrals

  • mo_num

  • n_states

  • one_e_dm_mo_alpha_for_dft

  • one_e_dm_mo_beta_for_dft

short_range_hartree

File : dft_utils_one_e/sr_coulomb.irp.f

double precision, allocatable   :: short_range_hartree_operator (mo_num,mo_num,N_states)
double precision, allocatable   :: short_range_hartree  (N_states)

short_range_Hartree_operator(i,j) = \(\int dr i(r)j(r) \int r' \rho(r') W_{ee}^{sr}\)

short_range_Hartree = \(1/2 \sum_{i,j} \rho_{ij} \mathtt{short_range_Hartree_operator}(i,j)\)

= \(1/2 \int dr \int r' \rho(r) \rho(r') W_{ee}^{sr}\)

Needs:

  • mo_integrals_erf_map

  • mo_integrals_map

  • mo_num

  • mo_two_e_integrals_erf_in_map

  • mo_two_e_integrals_in_map

  • n_states

  • one_e_dm_average_mo_for_dft

  • one_e_dm_mo_for_dft

Needed by:

  • effective_one_e_potential

  • trace_v_xc

short_range_hartree_operator

File : dft_utils_one_e/sr_coulomb.irp.f

double precision, allocatable   :: short_range_hartree_operator (mo_num,mo_num,N_states)
double precision, allocatable   :: short_range_hartree  (N_states)

short_range_Hartree_operator(i,j) = \(\int dr i(r)j(r) \int r' \rho(r') W_{ee}^{sr}\)

short_range_Hartree = \(1/2 \sum_{i,j} \rho_{ij} \mathtt{short_range_Hartree_operator}(i,j)\)

= \(1/2 \int dr \int r' \rho(r) \rho(r') W_{ee}^{sr}\)

Needs:

  • mo_integrals_erf_map

  • mo_integrals_map

  • mo_num

  • mo_two_e_integrals_erf_in_map

  • mo_two_e_integrals_in_map

  • n_states

  • one_e_dm_average_mo_for_dft

  • one_e_dm_mo_for_dft

Needed by:

  • effective_one_e_potential

  • trace_v_xc

Subroutines / functions

berf:

File : dft_utils_one_e/exc_sr_lda.irp.f

function berf(a)
dberfda:

File : dft_utils_one_e/exc_sr_lda.irp.f

function dberfda(a)
dpol:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function dpol(rs)
dpold:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function dpold(rs)
dpoldd:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function dpoldd(rs)
ec_lda:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ec_lda(rho_a,rho_b,ec,vc_a,vc_b)

Called by:

  • ec_pbe_only()

  • ec_pbe_sr()

  • energy_c_lda

Calls:

  • ecpw()

ec_lda_sr:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ec_lda_sr(mu,rho_a,rho_b,ec,vc_a,vc_b)

Called by:

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vxc_alpha_lda_w

  • ec_pbe_only()

  • ec_pbe_sr()

  • energy_c_sr_lda

  • energy_sr_x_lda

Calls:

  • ecorrlr()

  • ecpw()

  • vcorrlr()

ec_only_lda_sr:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ec_only_lda_sr(mu,rho_a,rho_b,ec)

Calls:

  • ecorrlr()

  • ecpw()

ec_pbe_only:

File : dft_utils_one_e/exc_sr_pbe.irp.f

subroutine ec_pbe_only(mu,rhoc,rhoo,sigmacc,sigmaco,sigmaoo,ec)

Short-range pbe correlation energy functional for erf interaction

input : ==========

mu = range separated parameter

rhoc, rhoo = total density and spin density

sigmacc = square of the gradient of the total density

sigmaco = square of the gradient of the spin density

sigmaoo = scalar product between the gradient of the total density and the one of the spin density

output: ==========

ec = correlation energy

Calls:

  • ec_lda()

  • ec_lda_sr()

ec_pbe_sr:

File : dft_utils_one_e/exc_sr_pbe.irp.f

subroutine ec_pbe_sr(mu,rhoc,rhoo,sigmacc,sigmaco,sigmaoo,ec,vrhoc,vrhoo,vsigmacc,vsigmaco,vsigmaoo)

Short-range pbe correlation energy functional for erf interaction

input : ==========

mu = range separated parameter

rhoc, rhoo = total density and spin density

sigmacc = square of the gradient of the total density

sigmaco = square of the gradient of the spin density

sigmaoo = scalar product between the gradient of the total density and the one of the spin density

output: ==========

ec = correlation energy

all variables v** are energy derivatives with respect to components of the density

vrhoc = derivative with respect to the total density

vrhoo = derivative with respect to spin density

vsigmacc = derivative with respect to the square of the gradient of the total density

vsigmaco = derivative with respect to scalar product between the gradients of total and spin densities

vsigmaoo = derivative with respect to the square of the gradient of the psin density

Called by:

  • gga_sr_type_functionals()

  • gga_type_functionals()

Calls:

  • ec_lda()

  • ec_lda_sr()

ecorrlr:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ecorrlr(rs,z,mu,eclr)

Called by:

  • ec_lda_sr()

  • ec_only_lda_sr()

Calls:

  • ecpw()

ecpw:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ecPW(x,y,ec,ecd,ecz,ecdd,eczd)

Called by:

  • ec_lda()

  • ec_lda_sr()

  • ec_only_lda_sr()

  • ecorrlr()

  • vcorrlr()

Calls:

  • gpw()

ex_lda:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ex_lda(rho_a,rho_b,ex,vx_a,vx_b)

Called by:

  • energy_x_lda

ex_lda_sr:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine ex_lda_sr(mu,rho_a,rho_b,ex,vx_a,vx_b)

Called by:

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vxc_alpha_lda_w

  • energy_sr_x_lda

  • energy_x_sr_lda

  • ex_pbe_sr()

  • ex_pbe_sr_only()

ex_pbe_sr:

File : dft_utils_one_e/exc_sr_pbe.irp.f

subroutine ex_pbe_sr(mu,rho_a,rho_b,grd_rho_a_2,grd_rho_b_2,grd_rho_a_b,ex,vx_rho_a,vx_rho_b,vx_grd_rho_a_2,vx_grd_rho_b_2,vx_grd_rho_a_b)

mu = range separation parameter rho_a = density alpha rho_b = density beta grd_rho_a_2 = (gradient rho_a)^2 grd_rho_b_2 = (gradient rho_b)^2 grd_rho_a_b = (gradient rho_a).(gradient rho_b) ex = exchange energy density at the density and corresponding gradients of the density vx_rho_a = d ex / d rho_a vx_rho_b = d ex / d rho_b vx_grd_rho_a_2 = d ex / d grd_rho_a_2 vx_grd_rho_b_2 = d ex / d grd_rho_b_2 vx_grd_rho_a_b = d ex / d grd_rho_a_b

Called by:

  • gga_sr_type_functionals()

  • gga_type_functionals()

Calls:

  • ex_lda_sr()

ex_pbe_sr_only:

File : dft_utils_one_e/exc_sr_pbe.irp.f

subroutine ex_pbe_sr_only(mu,rho_a,rho_b,grd_rho_a_2,grd_rho_b_2,grd_rho_a_b,ex)

rho_a = density alpha rho_b = density beta grd_rho_a_2 = (gradient rho_a)^2 grd_rho_b_2 = (gradient rho_b)^2 grd_rho_a_b = (gradient rho_a).(gradient rho_b) ex = exchange energy density at point r

Calls:

  • ex_lda_sr()

g0d:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function g0d(rs)
g0dd:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function g0dd(rs)
g0f:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function g0f(x)
gpw:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine GPW(x,Ac,alfa1,beta1,beta2,beta3,beta4,G,Gd,Gdd)

Called by:

  • ecpw()

grad_rho_ab_to_grad_rho_oc:

File : dft_utils_one_e/rho_ab_to_rho_tot.irp.f

subroutine grad_rho_ab_to_grad_rho_oc(grad_rho_a_2,grad_rho_b_2,grad_rho_a_b,grad_rho_o_2,grad_rho_c_2,grad_rho_o_c)

Called by:

  • gga_sr_type_functionals()

  • gga_type_functionals()

qrpa:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function Qrpa(x)
qrpad:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function Qrpad(x)
qrpadd:

File : dft_utils_one_e/exc_sr_lda.irp.f

double precision function Qrpadd(x)
rho_ab_to_rho_oc:

File : dft_utils_one_e/rho_ab_to_rho_tot.irp.f

subroutine rho_ab_to_rho_oc(rho_a,rho_b,rho_o,rho_c)

Called by:

  • gga_sr_type_functionals()

  • gga_type_functionals()

rho_oc_to_rho_ab:

File : dft_utils_one_e/rho_ab_to_rho_tot.irp.f

subroutine rho_oc_to_rho_ab(rho_o,rho_c,rho_a,rho_b)
v_grad_rho_oc_to_v_grad_rho_ab:

File : dft_utils_one_e/rho_ab_to_rho_tot.irp.f

subroutine v_grad_rho_oc_to_v_grad_rho_ab(v_grad_rho_o_2,v_grad_rho_c_2,v_grad_rho_o_c,v_grad_rho_a_2,v_grad_rho_b_2,v_grad_rho_a_b)

Called by:

  • gga_sr_type_functionals()

  • gga_type_functionals()

v_rho_ab_to_v_rho_oc:

File : dft_utils_one_e/rho_ab_to_rho_tot.irp.f

subroutine v_rho_ab_to_v_rho_oc(v_rho_a,v_rho_b,v_rho_o,v_rho_c)
v_rho_oc_to_v_rho_ab:

File : dft_utils_one_e/rho_ab_to_rho_tot.irp.f

subroutine v_rho_oc_to_v_rho_ab(v_rho_o,v_rho_c,v_rho_a,v_rho_b)

Called by:

  • gga_sr_type_functionals()

  • gga_type_functionals()

vcorrlr:

File : dft_utils_one_e/exc_sr_lda.irp.f

subroutine vcorrlr(rs,z,mu,vclrup,vclrdown,vclrupd,vclrdownd)

Called by:

  • ec_lda_sr()

Calls:

  • ecpw()