9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-11-01 19:23:38 +01:00
qp2/src/dft_utils_in_r/ao_in_r.irp.f
Anthony Scemama 49e9488f62
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
2019-02-22 19:19:58 +01:00

124 lines
3.8 KiB
Fortran

BEGIN_PROVIDER[double precision, aos_in_r_array, (ao_num,n_points_final_grid)]
&BEGIN_PROVIDER[double precision, aos_in_r_array_transp, (n_points_final_grid,ao_num)]
implicit none
BEGIN_DOC
! aos_in_r_array(i,j) = value of the ith ao on the jth grid point
!
! aos_in_r_array_transp(i,j) = value of the jth ao on the ith grid point
END_DOC
integer :: i,j
double precision :: aos_array(ao_num), r(3)
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)
call give_all_aos_at_r(r,aos_array)
do j = 1, ao_num
aos_in_r_array(j,i) = aos_array(j)
aos_in_r_array_transp(i,j) = aos_array(j)
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, aos_grad_in_r_array, (ao_num,n_points_final_grid,3)]
implicit none
BEGIN_DOC
! aos_grad_in_r_array(i,j,k) = value of the kth component of the gradient of ith ao on the jth grid point
!
! k = 1 : x, k= 2, y, k 3, z
END_DOC
integer :: i,j,m
double precision :: aos_array(ao_num), r(3)
double precision :: aos_grad_array(3,ao_num)
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)
call give_all_aos_and_grad_at_r(r,aos_array,aos_grad_array)
do m = 1, 3
do j = 1, ao_num
aos_grad_in_r_array(j,i,m) = aos_grad_array(m,j)
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, aos_grad_in_r_array_transp, (n_points_final_grid,ao_num,3)]
implicit none
BEGIN_DOC
! aos_grad_in_r_array_transp(i,j,k) = value of the kth component of the gradient of jth ao on the ith grid point
!
! k = 1 : x, k= 2, y, k 3, z
END_DOC
integer :: i,j,m
double precision :: aos_array(ao_num), r(3)
double precision :: aos_grad_array(3,ao_num)
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)
call give_all_aos_and_grad_at_r(r,aos_array,aos_grad_array)
do m = 1, 3
do j = 1, ao_num
aos_grad_in_r_array_transp(i,j,m) = aos_grad_array(m,j)
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, aos_grad_in_r_array_transp_xyz, (3,ao_num,n_points_final_grid)]
implicit none
BEGIN_DOC
! aos_grad_in_r_array_transp_xyz(k,i,j) = value of the kth component of the gradient of jth ao on the ith grid point
!
! k = 1 : x, k= 2, y, k 3, z
END_DOC
integer :: i,j,m
double precision :: aos_array(ao_num), r(3)
double precision :: aos_grad_array(3,ao_num)
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)
call give_all_aos_and_grad_at_r(r,aos_array,aos_grad_array)
do m = 1, 3
do j = 1, ao_num
aos_grad_in_r_array_transp_xyz(m,j,i) = aos_grad_array(m,j)
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER[double precision, aos_lapl_in_r_array, (ao_num,n_points_final_grid,3)]
&BEGIN_PROVIDER[double precision, aos_lapl_in_r_array_transp, (n_points_final_grid,ao_num,3)]
implicit none
BEGIN_DOC
! aos_lapl_in_r_array(i,j,k) = value of the kth component of the laplacian of ith ao on the jth grid point
!
! aos_lapl_in_r_array_transp(i,j,k) = value of the kth component of the laplacian of jth ao on the ith grid point
!
! k = 1 : x, k= 2, y, k 3, z
END_DOC
integer :: i,j,m
double precision :: aos_array(ao_num), r(3)
double precision :: aos_grad_array(ao_num,3)
double precision :: aos_lapl_array(ao_num,3)
do m = 1, 3
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)
call give_all_aos_and_grad_and_lapl_at_r(r,aos_array,aos_grad_array,aos_lapl_array)
do j = 1, ao_num
aos_lapl_in_r_array(j,i,m) = aos_lapl_array(j,m)
aos_lapl_in_r_array_transp(i,j,m) = aos_lapl_array(j,m)
enddo
enddo
enddo
END_PROVIDER