dft_utils_one_e¶
This module contains all the one-body related quantities needed to perform DFT or RS-DFT calculations. Therefore, it contains most of the properties which depends on the one-body density and density matrix.
The most important files and variables are:
- The general providers for the x/c energies in
e_xc_general.irp.f
- The general providers for the x/c potentials in
pot_general.irp.f
- The short-range hartree operator and all related quantities in
sr_coulomb.irp.f
These providers will be used in many DFT-related programs, such as ks_scf.irp.f
or rs_ks_scf.irp.f
.
It is also needed to compute the effective one-body operator needed in multi-determinant RS-DFT (see plugins by eginer).
Some other interesting quantities:
- The LDA and PBE providers for the x/c energies in
e_xc.irp.f
andsr_exc.irp.f
- The LDA and PBE providers for the x/c potentials on the AO basis in
pot_ao.irp.f
andsr_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
andexc_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
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
n_states
-
aos_dsr_vc_alpha_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_dsr_vc_beta_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_dsr_vx_alpha_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_dsr_vx_beta_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_dvc_alpha_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_dvc_beta_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_dvx_alpha_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_dvx_beta_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_sr_vc_alpha_lda_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_sr_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (sr_v^x_alpha(r_j) + sr_v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
final_grid_points
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_sr_vc_alpha_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_sr_vc_beta_lda_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_sr_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (sr_v^x_alpha(r_j) + sr_v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
final_grid_points
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_sr_vc_beta_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_sr_vx_alpha_lda_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_sr_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (sr_v^x_alpha(r_j) + sr_v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
final_grid_points
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_sr_vx_alpha_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_sr_vx_beta_lda_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_sr_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_sr_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (sr_v^x_alpha(r_j) + sr_v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
final_grid_points
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_sr_vx_beta_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
aos_vc_alpha_lda_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_vc_alpha_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_vc_beta_lda_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_vc_beta_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_vx_alpha_lda_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_vx_alpha_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
aos_vx_beta_lda_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vc_beta_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_alpha_lda_w (n_points_final_grid,ao_num,N_states) double precision, allocatable :: aos_vx_beta_lda_w (n_points_final_grid,ao_num,N_states)
aos_vxc_alpha_LDA_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_in_r_array
n_states
one_e_dm_alpha_at_r
Needed by:
-
aos_vx_beta_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
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_num
n_states
Needed by:
-
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_num
n_states
Needed by:
-
energy_c
¶ File :
dft_utils_one_e/e_xc_general.irp.f
double precision, allocatable :: energy_x (N_states) double precision, allocatable :: energy_c (N_states)
correlation and exchange energies general providers.
Needs:
correlation_functional
energy_sr_x_lda
energy_sr_x_pbe
exchange_functional
n_states
Needed by:
-
energy_c_lda
¶ File :
dft_utils_one_e/e_xc.irp.f
double precision, allocatable :: energy_x_lda (N_states) double precision, allocatable :: energy_c_lda (N_states)
exchange/correlation energy with the short range LDA functional
Needs:
n_states
-
energy_c_pbe
¶ File :
dft_utils_one_e/e_xc.irp.f
double precision, allocatable :: energy_x_pbe (N_states) double precision, allocatable :: energy_c_pbe (N_states)
exchange/correlation energy with the short range PBE functional
Needs:
correlation_functional
exchange_functional
n_states
one_e_dm_and_grad_alpha_in_r
-
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:
n_points_final_grid
n_states
Needed by:
-
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:
correlation_functional
exchange_functional
final_grid_points
n_states
one_e_dm_and_grad_alpha_in_r
Needed by:
-
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:
n_points_final_grid
n_states
Needed by:
-
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:
correlation_functional
exchange_functional
final_grid_points
n_states
one_e_dm_and_grad_alpha_in_r
Needed by:
-
energy_x
¶ File :
dft_utils_one_e/e_xc_general.irp.f
double precision, allocatable :: energy_x (N_states) double precision, allocatable :: energy_c (N_states)
correlation and exchange energies general providers.
Needs:
correlation_functional
energy_sr_x_lda
energy_sr_x_pbe
exchange_functional
n_states
Needed by:
-
energy_x_lda
¶ File :
dft_utils_one_e/e_xc.irp.f
double precision, allocatable :: energy_x_lda (N_states) double precision, allocatable :: energy_c_lda (N_states)
exchange/correlation energy with the short range LDA functional
Needs:
n_states
-
energy_x_pbe
¶ File :
dft_utils_one_e/e_xc.irp.f
double precision, allocatable :: energy_x_pbe (N_states) double precision, allocatable :: energy_c_pbe (N_states)
exchange/correlation energy with the short range PBE functional
Needs:
correlation_functional
exchange_functional
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
exchange_functional
correlation_functional
n_states
Called by:
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
exchange_functional
correlation_functional
Called by:
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()
-
grad_aos_dsr_vc_alpha_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
grad_aos_dsr_vc_beta_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
grad_aos_dsr_vx_alpha_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
grad_aos_dsr_vx_beta_pbe_w
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: aos_sr_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_sr_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dsr_vx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
ao_num
aos_grad_in_r_array
aos_in_r_array
correlation_functional
exchange_functional
final_grid_points
mu_erf_dft
Needed by:
-
grad_aos_dvc_alpha_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
grad_aos_dvc_beta_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
grad_aos_dvx_alpha_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
grad_aos_dvx_beta_pbe_w
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: aos_vc_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vc_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_alpha_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_vx_beta_pbe_w (ao_num,n_points_final_grid,N_states) double precision, allocatable :: aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvc_beta_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_alpha_pbe_w (ao_num,n_points_final_grid,3,N_states) double precision, allocatable :: grad_aos_dvx_beta_pbe_w (ao_num,n_points_final_grid,3,N_states)
aos_vxc_alpha_PBE_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
Needs:
correlation_functional
exchange_functional
final_grid_points
Needed by:
-
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:
-
potential_c_alpha_ao
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao (ao_num,ao_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the AO basis
Needs:
ao_num
correlation_functional
exchange_functional
Needed by:
-
potential_c_alpha_ao_lda
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_lda (ao_num,ao_num,N_states)
short range exchange/correlation alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_c_alpha_ao_pbe
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_c_alpha_mo
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_x_beta_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_beta_mo (mo_num,mo_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the MO basis
Needs:
ao_num
mo_coef
mo_num
n_states
Needed by:
-
potential_c_beta_ao
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao (ao_num,ao_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the AO basis
Needs:
ao_num
correlation_functional
exchange_functional
Needed by:
-
potential_c_beta_ao_lda
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_lda (ao_num,ao_num,N_states)
short range exchange/correlation alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_c_beta_ao_pbe
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_c_beta_mo
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_x_beta_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_beta_mo (mo_num,mo_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the MO basis
Needs:
ao_num
mo_coef
mo_num
n_states
Needed by:
-
potential_sr_c_alpha_ao_lda
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_c_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_beta_ao_lda (ao_num,ao_num,N_states)
short range correlation alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_sr_c_alpha_ao_pbe
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_sr_c_beta_ao_lda
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_c_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_beta_ao_lda (ao_num,ao_num,N_states)
short range correlation alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_sr_c_beta_ao_pbe
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_sr_x_alpha_ao_lda
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_x_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_x_beta_ao_lda (ao_num,ao_num,N_states)
short range exchange alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_sr_x_alpha_ao_pbe
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_sr_x_beta_ao_lda
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_x_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_x_beta_ao_lda (ao_num,ao_num,N_states)
short range exchange alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_sr_x_beta_ao_pbe
¶ File :
dft_utils_one_e/sr_pot_ao.irp.f
double precision, allocatable :: potential_sr_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_sr_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_x_alpha_ao
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao (ao_num,ao_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the AO basis
Needs:
ao_num
correlation_functional
exchange_functional
Needed by:
-
potential_x_alpha_ao_lda
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_lda (ao_num,ao_num,N_states)
short range exchange/correlation alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_x_alpha_ao_pbe
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_x_alpha_mo
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_x_beta_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_beta_mo (mo_num,mo_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the MO basis
Needs:
ao_num
mo_coef
mo_num
n_states
Needed by:
-
potential_x_beta_ao
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao (ao_num,ao_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the AO basis
Needs:
ao_num
correlation_functional
exchange_functional
Needed by:
-
potential_x_beta_ao_lda
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_lda (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_lda (ao_num,ao_num,N_states)
short range exchange/correlation alpha/beta potentials with LDA functional on the AO basis
Needs:
ao_num
aos_in_r_array
n_states
Needed by:
-
potential_x_beta_ao_pbe
¶ File :
dft_utils_one_e/pot_ao.irp.f
double precision, allocatable :: potential_x_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_x_beta_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_alpha_ao_pbe (ao_num,ao_num,N_states) double precision, allocatable :: potential_c_beta_ao_pbe (ao_num,ao_num,N_states)
exchange/correlation alpha/beta potentials with the short range PBE functional on the AO basis
Needs:
ao_num
aos_grad_in_r_array
n_points_final_grid
n_states
Needed by:
-
potential_x_beta_mo
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: potential_x_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_x_beta_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_alpha_mo (mo_num,mo_num,N_states) double precision, allocatable :: potential_c_beta_mo (mo_num,mo_num,N_states)
general providers for the alpha/beta exchange/correlation potentials on the MO basis
Needs:
ao_num
mo_coef
mo_num
n_states
Needed by:
-
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
-
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
-
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
-
shifting_constant
¶ File :
dft_utils_one_e/shifted_potential.irp.f
double precision, allocatable :: shifting_constant (N_states)
shifting_constant = (E_{Hxc} - <Psi | V_{Hxc} | Psi>) / N_elec constant to add to the potential in order to obtain the variational energy as the eigenvalue of the effective long-range Hamiltonian (see original paper of Levy PRL 113, 113002 (2014), equation (17) )
Needs:
n_states
short_range_hartree_operator
-
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:
Needed by:
-
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:
Needed by:
-
trace_v_h
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: trace_v_xc (N_states) double precision, allocatable :: trace_v_h (N_states) double precision, allocatable :: trace_v_hxc (N_states)
Trace_v_xc = sum_{i,j} (rho_{ij}_alpha v^{xc}_{ij}^alpha + rho_{ij}_beta v^{xc}_{ij}^beta) Trace_v_Hxc = sum_{i,j} v^{H}_{ij} (rho_{ij}_alpha + rho_{ij}_beta) Trace_v_Hxc = sum_{i,j} rho_{ij} v^{Hxc}_{ij}
Needs:
mo_num
n_states
Needed by:
-
trace_v_hxc
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: trace_v_xc (N_states) double precision, allocatable :: trace_v_h (N_states) double precision, allocatable :: trace_v_hxc (N_states)
Trace_v_xc = sum_{i,j} (rho_{ij}_alpha v^{xc}_{ij}^alpha + rho_{ij}_beta v^{xc}_{ij}^beta) Trace_v_Hxc = sum_{i,j} v^{H}_{ij} (rho_{ij}_alpha + rho_{ij}_beta) Trace_v_Hxc = sum_{i,j} rho_{ij} v^{Hxc}_{ij}
Needs:
mo_num
n_states
Needed by:
-
trace_v_xc
¶ File :
dft_utils_one_e/pot_general.irp.f
double precision, allocatable :: trace_v_xc (N_states) double precision, allocatable :: trace_v_h (N_states) double precision, allocatable :: trace_v_hxc (N_states)
Trace_v_xc = sum_{i,j} (rho_{ij}_alpha v^{xc}_{ij}^alpha + rho_{ij}_beta v^{xc}_{ij}^beta) Trace_v_Hxc = sum_{i,j} v^{H}_{ij} (rho_{ij}_alpha + rho_{ij}_beta) Trace_v_Hxc = sum_{i,j} rho_{ij} v^{Hxc}_{ij}
Needs:
mo_num
n_states
Needed by:
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()
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:
ec_pbe_only()
ec_pbe_sr()
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:
-
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:
energy_sr_x_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()