9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-06-19 11:02:04 +02:00
qp2/src/functionals/sr_lda.irp.f
Anthony Scemama 8b22e38c9c
Develop (#15)
* fixed laplacian of aos

* corrected the laplacians of aos

* added dft_one_e

* added new feature for new dft functionals

* changed the configure to add new functionals

* changed the configure

* added dft_one_e/README.rst

* added README.rst in new_functionals

* added source/programmers_guide/new_ks.rst

* Thesis Yann

* Added gmp installation in configure

* improved qp_e_conv_fci

* Doc

* Typos

* Added variance_max

* Fixed completion in qp_create

* modif TODO

* fixed DFT potential for n_states gt 1

* improved pot pbe

* trying to improve sr PBE

* fixed potential pbe

* fixed the vxc smashed for pbe sr and normal

* Comments in selection

* bug fixed by peter

* Fixed bug with zero beta electrons

* Update README.rst

* Update e_xc_new_func.irp.f

* Update links.rst

* Update quickstart.rst

* Update quickstart.rst

* updated cipsi

* Fixed energies of non-expected s2 (#9)

* Moved diag_algorithm in Davdison

* Add print_ci_vector in tools (#11)

* Fixed energies of non-expected s2

* Moved diag_algorithm in Davdison

* Fixed travis

* Added print_ci_vector

* Documentation

* Cleaned qp_set_mo_class.ml

* Removed Core in taskserver

* Merge develop-toto and manus (#12)

* Fixed energies of non-expected s2

* Moved diag_algorithm in Davdison

* Fixed travis

* Added print_ci_vector

* Documentation

* Cleaned qp_set_mo_class.ml

* Removed Core in taskserver

* Frozen core for heavy atoms

* Improved molden module

* In sync with manus

* Fixed some of the documentation errors

* Develop toto (#13)

* Fixed energies of non-expected s2

* Moved diag_algorithm in Davdison

* Fixed travis

* Added print_ci_vector

* Documentation

* Cleaned qp_set_mo_class.ml

* Removed Core in taskserver

* Frozen core for heavy atoms

* Improved molden module

* In sync with manus

* Fixed some of the documentation errors

* Develop manus (#14)

* modified printing for rpt2

* Comment

* Fixed plugins

* Scripting for functionals

* Documentation

* Develop (#10)

* fixed laplacian of aos

* corrected the laplacians of aos

* added dft_one_e

* added new feature for new dft functionals

* changed the configure to add new functionals

* changed the configure

* added dft_one_e/README.rst

* added README.rst in new_functionals

* added source/programmers_guide/new_ks.rst

* Thesis Yann

* Added gmp installation in configure

* improved qp_e_conv_fci

* Doc

* Typos

* Added variance_max

* Fixed completion in qp_create

* modif TODO

* fixed DFT potential for n_states gt 1

* improved pot pbe

* trying to improve sr PBE

* fixed potential pbe

* fixed the vxc smashed for pbe sr and normal

* Comments in selection

* bug fixed by peter

* Fixed bug with zero beta electrons

* Update README.rst

* Update e_xc_new_func.irp.f

* Update links.rst

* Update quickstart.rst

* Update quickstart.rst

* updated cipsi

* Fixed energies of non-expected s2 (#9)

* Moved diag_algorithm in Davdison

* some modifs

* modified gfortran_debug.cfg

* fixed automatization of functionals

* modified e_xc_general.irp.f

* minor modifs in ref_bitmask.irp.f

* modifying functionals

* rs_ks_scf and ks_scf compiles with the automatic handling of functionals

* removed prints

* fixed configure

* fixed the new functionals

* Merge toto

* modified automatic functionals

* Changed python into python2

* from_xyz suppressed

* Cleaning repo

* Update README.md

* Update README.md

* Contributors

* Update GITHUB.md

* bibtex
2019-03-07 16:29:06 +01:00

196 lines
7.4 KiB
Fortran

BEGIN_PROVIDER[double precision, energy_x_sr_lda, (N_states) ]
implicit none
BEGIN_DOC
! exchange energy with the short range lda functional
END_DOC
integer :: istate,i,j
double precision :: r(3)
double precision :: mu,weight
double precision :: e_x,vx_a,vx_b
double precision, allocatable :: rhoa(:),rhob(:)
allocate(rhoa(N_states), rhob(N_states))
energy_x_sr_lda = 0.d0
do istate = 1, N_states
do i = 1, n_points_final_grid
r(1) = final_grid_points(1,i)
r(2) = final_grid_points(2,i)
r(3) = final_grid_points(3,i)
weight = final_weight_at_r_vector(i)
rhoa(istate) = one_e_dm_alpha_at_r(i,istate)
rhob(istate) = one_e_dm_beta_at_r(i,istate)
call ex_lda_sr(mu_erf_dft,rhoa(istate),rhob(istate),e_x,vx_a,vx_b)
energy_x_sr_lda(istate) += weight * e_x
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, energy_c_sr_lda, (N_states) ]
implicit none
BEGIN_DOC
! exchange energy with the short range lda functional
END_DOC
integer :: istate,i,j
double precision :: r(3)
double precision :: mu,weight
double precision :: e_c,vc_a,vc_b
double precision, allocatable :: rhoa(:),rhob(:)
allocate(rhoa(N_states), rhob(N_states))
energy_c_sr_lda = 0.d0
do istate = 1, N_states
do i = 1, n_points_final_grid
r(1) = final_grid_points(1,i)
r(2) = final_grid_points(2,i)
r(3) = final_grid_points(3,i)
weight = final_weight_at_r_vector(i)
rhoa(istate) = one_e_dm_alpha_at_r(i,istate)
rhob(istate) = one_e_dm_beta_at_r(i,istate)
call ec_lda_sr(mu_erf_dft,rhoa(istate),rhob(istate),e_c,vc_a,vc_b)
energy_c_sr_lda(istate) += weight * e_c
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, potential_x_alpha_ao_sr_lda,(ao_num,ao_num,N_states)]
&BEGIN_PROVIDER [double precision, potential_x_beta_ao_sr_lda,(ao_num,ao_num,N_states)]
implicit none
BEGIN_DOC
! short range exchange alpha/beta potentials with lda functional on the |AO| basis
END_DOC
! Second dimension is given as ao_num * N_states so that Lapack does the loop over N_states.
integer :: istate
do istate = 1, N_states
call dgemm('N','T',ao_num,ao_num,n_points_final_grid,1.d0, &
aos_in_r_array,size(aos_in_r_array,1), &
aos_sr_vx_alpha_lda_w,size(aos_sr_vx_alpha_lda_w,1),0.d0,&
potential_x_alpha_ao_sr_lda,size(potential_x_alpha_ao_sr_lda,1))
call dgemm('N','T',ao_num,ao_num,n_points_final_grid,1.d0, &
aos_in_r_array,size(aos_in_r_array,1), &
aos_sr_vx_beta_lda_w(1,1,istate),size(aos_sr_vx_beta_lda_w,1),0.d0,&
potential_x_beta_ao_sr_lda(1,1,istate),size(potential_x_beta_ao_sr_lda,1))
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, potential_c_alpha_ao_sr_lda,(ao_num,ao_num,N_states)]
&BEGIN_PROVIDER [double precision, potential_c_beta_ao_sr_lda,(ao_num,ao_num,N_states)]
implicit none
BEGIN_DOC
! short range correlation alpha/beta potentials with lda functional on the |AO| basis
END_DOC
! Second dimension is given as ao_num * N_states so that Lapack does the loop over N_states.
integer :: istate
do istate = 1, N_states
call dgemm('N','T',ao_num,ao_num,n_points_final_grid,1.d0, &
aos_in_r_array,size(aos_in_r_array,1), &
aos_sr_vc_alpha_lda_w(1,1,istate),size(aos_sr_vc_alpha_lda_w,1),0.d0,&
potential_c_alpha_ao_sr_lda(1,1,istate),size(potential_c_alpha_ao_sr_lda,1))
call dgemm('N','T',ao_num,ao_num,n_points_final_grid,1.d0, &
aos_in_r_array,size(aos_in_r_array,1), &
aos_sr_vc_beta_lda_w(1,1,istate),size(aos_sr_vc_beta_lda_w,1),0.d0,&
potential_c_beta_ao_sr_lda(1,1,istate),size(potential_c_beta_ao_sr_lda,1))
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, aos_sr_vc_alpha_lda_w, (ao_num,n_points_final_grid,N_states)]
&BEGIN_PROVIDER[double precision, aos_sr_vc_beta_lda_w, (ao_num,n_points_final_grid,N_states)]
&BEGIN_PROVIDER[double precision, aos_sr_vx_alpha_lda_w, (ao_num,n_points_final_grid,N_states)]
&BEGIN_PROVIDER[double precision, aos_sr_vx_beta_lda_w, (ao_num,n_points_final_grid,N_states)]
implicit none
BEGIN_DOC
! 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)
END_DOC
integer :: istate,i,j
double precision :: r(3)
double precision :: mu,weight
double precision :: e_c,sr_vc_a,sr_vc_b,e_x,sr_vx_a,sr_vx_b
double precision, allocatable :: rhoa(:),rhob(:)
allocate(rhoa(N_states), rhob(N_states))
do istate = 1, N_states
do i = 1, n_points_final_grid
r(1) = final_grid_points(1,i)
r(2) = final_grid_points(2,i)
r(3) = final_grid_points(3,i)
weight = final_weight_at_r_vector(i)
rhoa(istate) = one_e_dm_alpha_at_r(i,istate)
rhob(istate) = one_e_dm_beta_at_r(i,istate)
call ec_lda_sr(mu_erf_dft,rhoa(istate),rhob(istate),e_c,sr_vc_a,sr_vc_b)
call ex_lda_sr(mu_erf_dft,rhoa(istate),rhob(istate),e_x,sr_vx_a,sr_vx_b)
do j =1, ao_num
aos_sr_vc_alpha_lda_w(j,i,istate) = sr_vc_a * aos_in_r_array(j,i)*weight
aos_sr_vc_beta_lda_w(j,i,istate) = sr_vc_b * aos_in_r_array(j,i)*weight
aos_sr_vx_alpha_lda_w(j,i,istate) = sr_vx_a * aos_in_r_array(j,i)*weight
aos_sr_vx_beta_lda_w(j,i,istate) = sr_vx_b * aos_in_r_array(j,i)*weight
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, aos_sr_vxc_alpha_lda_w, (ao_num,n_points_final_grid,N_states)]
&BEGIN_PROVIDER[double precision, aos_sr_vxc_beta_lda_w, (ao_num,n_points_final_grid,N_states)]
implicit none
BEGIN_DOC
! aos_sr_vxc_alpha_lda_w(j,i) = ao_i(r_j) * (v^x_alpha(r_j) + v^c_alpha(r_j)) * W(r_j)
END_DOC
integer :: istate,i,j
double precision :: r(3)
double precision :: mu,weight
double precision :: e_c,sr_vc_a,sr_vc_b,e_x,sr_vx_a,sr_vx_b
double precision, allocatable :: rhoa(:),rhob(:)
double precision :: mu_local
mu_local = mu_erf_dft
allocate(rhoa(N_states), rhob(N_states))
do istate = 1, N_states
do i = 1, n_points_final_grid
r(1) = final_grid_points(1,i)
r(2) = final_grid_points(2,i)
r(3) = final_grid_points(3,i)
weight = final_weight_at_r_vector(i)
rhoa(istate) = one_e_dm_alpha_at_r(i,istate)
rhob(istate) = one_e_dm_beta_at_r(i,istate)
call ec_lda_sr(mu_local,rhoa(istate),rhob(istate),e_c,sr_vc_a,sr_vc_b)
call ex_lda_sr(mu_local,rhoa(istate),rhob(istate),e_x,sr_vx_a,sr_vx_b)
do j =1, ao_num
aos_sr_vxc_alpha_lda_w(j,i,istate) = (sr_vc_a + sr_vx_a) * aos_in_r_array(j,i)*weight
aos_sr_vxc_beta_lda_w(j,i,istate) = (sr_vc_b + sr_vx_b) * aos_in_r_array(j,i)*weight
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, potential_xc_alpha_ao_sr_lda,(ao_num,ao_num,N_states)]
&BEGIN_PROVIDER [double precision, potential_xc_beta_ao_sr_lda ,(ao_num,ao_num,N_states)]
implicit none
BEGIN_DOC
! short range exchange/correlation alpha/beta potentials with lda functional on the AO basis
END_DOC
integer :: istate
do istate = 1, N_states
call dgemm('N','T',ao_num,ao_num,n_points_final_grid,1.d0, &
aos_in_r_array,size(aos_in_r_array,1), &
aos_sr_vxc_alpha_lda_w(1,1,istate),size(aos_sr_vxc_alpha_lda_w,1),0.d0,&
potential_xc_alpha_ao_sr_lda(1,1,istate),size(potential_xc_alpha_ao_sr_lda,1))
call dgemm('N','T',ao_num,ao_num,n_points_final_grid,1.d0, &
aos_in_r_array,size(aos_in_r_array,1), &
aos_sr_vxc_beta_lda_w(1,1,istate),size(aos_sr_vxc_beta_lda_w,1),0.d0,&
potential_xc_beta_ao_sr_lda(1,1,istate),size(potential_xc_beta_ao_sr_lda,1))
enddo
END_PROVIDER