2019-07-23 15:12:25 +02:00
|
|
|
|
|
|
|
BEGIN_PROVIDER [ double precision, ao_pseudo_grid, (ao_num,-pseudo_lmax:pseudo_lmax,0:pseudo_lmax,nucl_num,pseudo_grid_size) ]
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! Grid points for f(|r-r_A|) = \int Y_{lm}^{C} (|r-r_C|, \Omega_C) \chi_i^{A} (r-r_A) d\Omega_C
|
|
|
|
!
|
|
|
|
! <img src="http://latex.codecogs.com/gif.latex?f(|r-r_A|)&space;=&space;\int&space;Y_{lm}^{C}&space;(|r-r_C|,&space;\Omega_C)&space;\chi_i^{A}&space;(r-r_A)&space;d\Omega_C"
|
|
|
|
! title="f(|r-r_A|) = \int Y_{lm}^{C} (|r-r_C|, \Omega_C) \chi_i^{A} (r-r_A) d\Omega_C" />
|
|
|
|
END_DOC
|
|
|
|
! l,m : Y(l,m) parameters
|
|
|
|
! c(3) : pseudopotential center
|
|
|
|
! a(3) : Atomic Orbital center
|
|
|
|
! n_a(3) : Powers of x,y,z in the Atomic Orbital
|
|
|
|
! g_a : Atomic Orbital exponent
|
|
|
|
! r : Distance between the Atomic Orbital center and the considered point
|
|
|
|
double precision, external :: ylm_orb
|
|
|
|
integer :: n_a(3)
|
|
|
|
double precision :: a(3), c(3), g_a
|
|
|
|
integer :: i,j,k,l,m,n,p
|
|
|
|
double precision :: dr, Ulc
|
|
|
|
double precision :: y
|
|
|
|
double precision, allocatable :: r(:)
|
|
|
|
|
|
|
|
allocate (r(pseudo_grid_size))
|
|
|
|
dr = pseudo_grid_rmax/dble(pseudo_grid_size)
|
|
|
|
r(1) = 0.d0
|
|
|
|
do j=2,pseudo_grid_size
|
|
|
|
r(j) = r(j-1) + dr
|
|
|
|
enddo
|
|
|
|
|
|
|
|
ao_pseudo_grid = 0.d0
|
2023-11-10 14:23:26 +01:00
|
|
|
!$OMP PARALLEL DO DEFAULT(SHARED) &
|
|
|
|
!$OMP PRIVATE(j,k,c,l,i,a,n_a,g_a,m,y)
|
2019-07-23 15:12:25 +02:00
|
|
|
do j=1,pseudo_grid_size
|
|
|
|
do k=1,nucl_num
|
|
|
|
c(1:3) = nucl_coord(k,1:3)
|
|
|
|
do l=0,pseudo_lmax
|
|
|
|
do i=1,ao_num
|
|
|
|
a(1:3) = nucl_coord(ao_nucl(i),1:3)
|
|
|
|
n_a(1:3) = ao_power(i,1:3)
|
|
|
|
do p=1,ao_prim_num(i)
|
|
|
|
g_a = ao_expo_ordered_transp(p,i)
|
|
|
|
do m=-l,l
|
|
|
|
y = ylm_orb(l,m,c,a,n_a,g_a,r(j))
|
|
|
|
ao_pseudo_grid(i,m,l,k,j) = ao_pseudo_grid(i,m,l,k,j) + &
|
|
|
|
ao_coef_normalized_ordered_transp(p,i)*y
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2023-11-10 14:23:26 +01:00
|
|
|
!$OMP END PARALLEL DO
|
2019-07-23 15:12:25 +02:00
|
|
|
deallocate(r)
|
|
|
|
|
|
|
|
END_PROVIDER
|
|
|
|
|
|
|
|
|
|
|
|
BEGIN_PROVIDER [ double precision, mo_pseudo_grid, (ao_num,-pseudo_lmax:pseudo_lmax,0:pseudo_lmax,nucl_num,pseudo_grid_size) ]
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! Grid points for f(|r-r_A|) = \int Y_{lm}^{C} (|r-r_C|, \Omega_C) \phi_i^{A} (r-r_A) d\Omega_C
|
|
|
|
!
|
|
|
|
! <img src="http://latex.codecogs.com/gif.latex?f(|r-r_A|)&space;=&space;\int&space;Y_{lm}^{C}&space;(|r-r_C|,&space;\Omega_C)&space;\chi_i^{A}&space;(r-r_A)&space;d\Omega_C"
|
|
|
|
! title="f(|r-r_A|) = \int Y_{lm}^{C} (|r-r_C|, \Omega_C) \chi_i^{A} (r-r_A) d\Omega_C" />
|
|
|
|
END_DOC
|
|
|
|
! l,m : Y(l,m) parameters
|
|
|
|
! c(3) : pseudopotential center
|
|
|
|
! a(3) : Atomic Orbital center
|
|
|
|
! n_a(3) : Powers of x,y,z in the Atomic Orbital
|
|
|
|
! g_a : Atomic Orbital exponent
|
|
|
|
! r : Distance between the Atomic Orbital center and the considered point
|
|
|
|
double precision, external :: ylm_orb
|
|
|
|
integer :: n_a(3)
|
|
|
|
double precision :: a(3), c(3), g_a
|
|
|
|
integer :: i,j,k,l,m,n,p
|
|
|
|
double precision :: dr, Ulc
|
|
|
|
double precision :: y
|
|
|
|
double precision, allocatable :: r(:)
|
|
|
|
|
|
|
|
|
2023-11-10 14:23:26 +01:00
|
|
|
mo_pseudo_grid = 0.d0
|
|
|
|
call dgemm('T','N', mo_num, (2*pseudo_lmax+1)*(pseudo_lmax+1)*nucl_num*pseudo_grid_size, &
|
|
|
|
ao_num, 1.d0, mo_coef, size(mo_coef,1), ao_pseudo_grid, &
|
|
|
|
size(ao_pseudo_grid,1), 0.d0, mo_pseudo_grid, size(mo_pseudo_grid,1))
|
2019-07-23 15:12:25 +02:00
|
|
|
|
|
|
|
END_PROVIDER
|
|
|
|
|