9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-10-16 04:31:32 +02:00

cos x GTOs integ added

This commit is contained in:
Abdallah Ammar 2023-03-04 17:49:48 +01:00
parent 4a23b63a42
commit 3c161384cc
22 changed files with 4591 additions and 674 deletions

View File

@ -12,21 +12,21 @@ double precision function ao_value(i,r)
integer :: power_ao(3)
double precision :: accu,dx,dy,dz,r2
num_ao = ao_nucl(i)
! power_ao(1:3)= ao_power(i,1:3)
! center_ao(1:3) = nucl_coord(num_ao,1:3)
! dx = (r(1) - center_ao(1))
! dy = (r(2) - center_ao(2))
! dz = (r(3) - center_ao(3))
! r2 = dx*dx + dy*dy + dz*dz
! dx = dx**power_ao(1)
! dy = dy**power_ao(2)
! dz = dz**power_ao(3)
power_ao(1:3)= ao_power(i,1:3)
center_ao(1:3) = nucl_coord(num_ao,1:3)
dx = (r(1) - center_ao(1))
dy = (r(2) - center_ao(2))
dz = (r(3) - center_ao(3))
r2 = dx*dx + dy*dy + dz*dz
dx = dx**power_ao(1)
dy = dy**power_ao(2)
dz = dz**power_ao(3)
accu = 0.d0
! do m=1,ao_prim_num(i)
! beta = ao_expo_ordered_transp(m,i)
! accu += ao_coef_normalized_ordered_transp(m,i) * dexp(-beta*r2)
! enddo
do m=1,ao_prim_num(i)
beta = ao_expo_ordered_transp(m,i)
accu += ao_coef_normalized_ordered_transp(m,i) * dexp(-beta*r2)
enddo
ao_value = accu * dx * dy * dz
end

View File

@ -1,2 +1,3 @@
ao_basis
pseudo
cosgtos_ao_int

View File

@ -1,75 +1,99 @@
BEGIN_PROVIDER [ double precision, ao_overlap,(ao_num,ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_overlap_x,(ao_num,ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_overlap_y,(ao_num,ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_overlap_z,(ao_num,ao_num) ]
implicit none
! ---
BEGIN_PROVIDER [ double precision, ao_overlap , (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_overlap_x, (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_overlap_y, (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_overlap_z, (ao_num, ao_num) ]
BEGIN_DOC
! Overlap between atomic basis functions:
!
! :math:`\int \chi_i(r) \chi_j(r) dr`
! Overlap between atomic basis functions:
!
! :math:`\int \chi_i(r) \chi_j(r) dr`
END_DOC
integer :: i,j,n,l
double precision :: f
integer :: dim1
implicit none
integer :: i, j, n, l, dim1, power_A(3), power_B(3)
double precision :: overlap, overlap_x, overlap_y, overlap_z
double precision :: alpha, beta, c
double precision :: A_center(3), B_center(3)
integer :: power_A(3), power_B(3)
ao_overlap = 0.d0
ao_overlap = 0.d0
ao_overlap_x = 0.d0
ao_overlap_y = 0.d0
ao_overlap_z = 0.d0
if (read_ao_integrals_overlap) then
call ezfio_get_ao_one_e_ints_ao_integrals_overlap(ao_overlap(1:ao_num, 1:ao_num))
print *, 'AO overlap integrals read from disk'
if(read_ao_integrals_overlap) then
call ezfio_get_ao_one_e_ints_ao_integrals_overlap(ao_overlap(1:ao_num, 1:ao_num))
print *, 'AO overlap integrals read from disk'
else
dim1=100
!$OMP PARALLEL DO SCHEDULE(GUIDED) &
!$OMP DEFAULT(NONE) &
!$OMP PRIVATE(A_center,B_center,power_A,power_B,&
!$OMP overlap_x,overlap_y, overlap_z, overlap, &
!$OMP alpha, beta,i,j,c) &
!$OMP SHARED(nucl_coord,ao_power,ao_prim_num, &
!$OMP ao_overlap_x,ao_overlap_y,ao_overlap_z,ao_overlap,ao_num,ao_coef_normalized_ordered_transp,ao_nucl, &
!$OMP ao_expo_ordered_transp,dim1)
do j=1,ao_num
A_center(1) = nucl_coord( ao_nucl(j), 1 )
A_center(2) = nucl_coord( ao_nucl(j), 2 )
A_center(3) = nucl_coord( ao_nucl(j), 3 )
power_A(1) = ao_power( j, 1 )
power_A(2) = ao_power( j, 2 )
power_A(3) = ao_power( j, 3 )
do i= 1,ao_num
B_center(1) = nucl_coord( ao_nucl(i), 1 )
B_center(2) = nucl_coord( ao_nucl(i), 2 )
B_center(3) = nucl_coord( ao_nucl(i), 3 )
power_B(1) = ao_power( i, 1 )
power_B(2) = ao_power( i, 2 )
power_B(3) = ao_power( i, 3 )
do n = 1,ao_prim_num(j)
alpha = ao_expo_ordered_transp(n,j)
do l = 1, ao_prim_num(i)
beta = ao_expo_ordered_transp(l,i)
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_x,overlap_y,overlap_z,overlap,dim1)
c = ao_coef_normalized_ordered_transp(n,j) * ao_coef_normalized_ordered_transp(l,i)
ao_overlap(i,j) += c * overlap
if(isnan(ao_overlap(i,j)))then
print*,'i,j',i,j
print*,'l,n',l,n
print*,'c,overlap',c,overlap
print*,overlap_x,overlap_y,overlap_z
stop
endif
ao_overlap_x(i,j) += c * overlap_x
ao_overlap_y(i,j) += c * overlap_y
ao_overlap_z(i,j) += c * overlap_z
if(use_cosgtos) then
!print*, ' use_cosgtos for ao_overlap ?', use_cosgtos
do j = 1, ao_num
do i = 1, ao_num
ao_overlap (i,j) = ao_overlap_cosgtos (i,j)
ao_overlap_x(i,j) = ao_overlap_cosgtos_x(i,j)
ao_overlap_y(i,j) = ao_overlap_cosgtos_y(i,j)
ao_overlap_z(i,j) = ao_overlap_cosgtos_z(i,j)
enddo
enddo
else
dim1=100
!$OMP PARALLEL DO SCHEDULE(GUIDED) &
!$OMP DEFAULT(NONE) &
!$OMP PRIVATE(A_center,B_center,power_A,power_B,&
!$OMP overlap_x,overlap_y, overlap_z, overlap, &
!$OMP alpha, beta,i,j,c) &
!$OMP SHARED(nucl_coord,ao_power,ao_prim_num, &
!$OMP ao_overlap_x,ao_overlap_y,ao_overlap_z,ao_overlap,ao_num,ao_coef_normalized_ordered_transp,ao_nucl, &
!$OMP ao_expo_ordered_transp,dim1)
do j=1,ao_num
A_center(1) = nucl_coord( ao_nucl(j), 1 )
A_center(2) = nucl_coord( ao_nucl(j), 2 )
A_center(3) = nucl_coord( ao_nucl(j), 3 )
power_A(1) = ao_power( j, 1 )
power_A(2) = ao_power( j, 2 )
power_A(3) = ao_power( j, 3 )
do i= 1,ao_num
B_center(1) = nucl_coord( ao_nucl(i), 1 )
B_center(2) = nucl_coord( ao_nucl(i), 2 )
B_center(3) = nucl_coord( ao_nucl(i), 3 )
power_B(1) = ao_power( i, 1 )
power_B(2) = ao_power( i, 2 )
power_B(3) = ao_power( i, 3 )
do n = 1,ao_prim_num(j)
alpha = ao_expo_ordered_transp(n,j)
do l = 1, ao_prim_num(i)
beta = ao_expo_ordered_transp(l,i)
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_x,overlap_y,overlap_z,overlap,dim1)
c = ao_coef_normalized_ordered_transp(n,j) * ao_coef_normalized_ordered_transp(l,i)
ao_overlap(i,j) += c * overlap
if(isnan(ao_overlap(i,j)))then
print*,'i,j',i,j
print*,'l,n',l,n
print*,'c,overlap',c,overlap
print*,overlap_x,overlap_y,overlap_z
stop
endif
ao_overlap_x(i,j) += c * overlap_x
ao_overlap_y(i,j) += c * overlap_y
ao_overlap_z(i,j) += c * overlap_z
enddo
enddo
enddo
enddo
enddo
enddo
!$OMP END PARALLEL DO
!$OMP END PARALLEL DO
endif
endif
if (write_ao_integrals_overlap) then
call ezfio_set_ao_one_e_ints_ao_integrals_overlap(ao_overlap(1:ao_num, 1:ao_num))
print *, 'AO overlap integrals written to disk'
@ -77,6 +101,8 @@
END_PROVIDER
! ---
BEGIN_PROVIDER [ double precision, ao_overlap_imag, (ao_num, ao_num) ]
implicit none
BEGIN_DOC
@ -85,6 +111,8 @@ BEGIN_PROVIDER [ double precision, ao_overlap_imag, (ao_num, ao_num) ]
ao_overlap_imag = 0.d0
END_PROVIDER
! ---
BEGIN_PROVIDER [ complex*16, ao_overlap_complex, (ao_num, ao_num) ]
implicit none
BEGIN_DOC
@ -98,41 +126,43 @@ BEGIN_PROVIDER [ complex*16, ao_overlap_complex, (ao_num, ao_num) ]
enddo
END_PROVIDER
! ---
BEGIN_PROVIDER [ double precision, ao_overlap_abs, (ao_num, ao_num) ]
BEGIN_PROVIDER [ double precision, ao_overlap_abs,(ao_num,ao_num) ]
implicit none
BEGIN_DOC
! Overlap between absolute values of atomic basis functions:
!
! :math:`\int |\chi_i(r)| |\chi_j(r)| dr`
! Overlap between absolute values of atomic basis functions:
!
! :math:`\int |\chi_i(r)| |\chi_j(r)| dr`
END_DOC
integer :: i,j,n,l
double precision :: f
integer :: dim1
double precision :: overlap, overlap_x, overlap_y, overlap_z
implicit none
integer :: i, j, n, l, dim1, power_A(3), power_B(3)
double precision :: overlap_x, overlap_y, overlap_z
double precision :: alpha, beta
double precision :: A_center(3), B_center(3)
integer :: power_A(3), power_B(3)
double precision :: lower_exp_val, dx
if (is_periodic) then
do j=1,ao_num
do i= 1,ao_num
ao_overlap_abs(i,j)= cdabs(ao_overlap_complex(i,j))
if(is_periodic) then
do j = 1, ao_num
do i = 1, ao_num
ao_overlap_abs(i,j) = cdabs(ao_overlap_complex(i,j))
enddo
enddo
else
dim1=100
lower_exp_val = 40.d0
!$OMP PARALLEL DO SCHEDULE(GUIDED) &
!$OMP DEFAULT(NONE) &
!$OMP PRIVATE(A_center,B_center,power_A,power_B, &
!$OMP overlap_x,overlap_y, overlap_z, overlap, &
!$OMP alpha, beta,i,j,dx) &
!$OMP SHARED(nucl_coord,ao_power,ao_prim_num, &
!$OMP ao_overlap_abs,ao_num,ao_coef_normalized_ordered_transp,ao_nucl,&
!$OMP ao_expo_ordered_transp,dim1,lower_exp_val)
!$OMP PARALLEL DO SCHEDULE(GUIDED) &
!$OMP DEFAULT(NONE) &
!$OMP PRIVATE(A_center,B_center,power_A,power_B, &
!$OMP overlap_x,overlap_y, overlap_z, &
!$OMP alpha, beta,i,j,dx) &
!$OMP SHARED(nucl_coord,ao_power,ao_prim_num, &
!$OMP ao_overlap_abs,ao_num,ao_coef_normalized_ordered_transp,ao_nucl,&
!$OMP ao_expo_ordered_transp,dim1,lower_exp_val)
do j=1,ao_num
A_center(1) = nucl_coord( ao_nucl(j), 1 )
A_center(2) = nucl_coord( ao_nucl(j), 2 )
@ -160,10 +190,14 @@ BEGIN_PROVIDER [ double precision, ao_overlap_abs,(ao_num,ao_num) ]
enddo
enddo
enddo
!$OMP END PARALLEL DO
!$OMP END PARALLEL DO
endif
END_PROVIDER
! ---
BEGIN_PROVIDER [ double precision, S_inv,(ao_num,ao_num) ]
implicit none
BEGIN_DOC

View File

@ -1,7 +1,10 @@
BEGIN_PROVIDER [ double precision, ao_deriv2_x,(ao_num,ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_deriv2_y,(ao_num,ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_deriv2_z,(ao_num,ao_num) ]
implicit none
! ---
BEGIN_PROVIDER [ double precision, ao_deriv2_x, (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_deriv2_y, (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_deriv2_z, (ao_num, ao_num) ]
BEGIN_DOC
! Second derivative matrix elements in the |AO| basis.
!
@ -11,114 +14,131 @@
! \langle \chi_i(x,y,z) | \frac{\partial^2}{\partial x^2} |\chi_j (x,y,z) \rangle
!
END_DOC
integer :: i,j,n,l
double precision :: f
integer :: dim1
implicit none
integer :: i, j, n, l, dim1, power_A(3), power_B(3)
double precision :: overlap, overlap_y, overlap_z
double precision :: overlap_x0, overlap_y0, overlap_z0
double precision :: alpha, beta, c
double precision :: A_center(3), B_center(3)
integer :: power_A(3), power_B(3)
double precision :: d_a_2,d_2
dim1=100
! -- Dummy call to provide everything
A_center(:) = 0.d0
B_center(:) = 1.d0
alpha = 1.d0
beta = .1d0
power_A = 1
power_B = 0
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,d_a_2,overlap_z,overlap,dim1)
! --
if(use_cosgtos) then
!print*, 'use_cosgtos for ao_kinetic_integrals ?', use_cosgtos
!$OMP PARALLEL DO SCHEDULE(GUIDED) &
!$OMP DEFAULT(NONE) &
!$OMP PRIVATE(A_center,B_center,power_A,power_B,&
!$OMP overlap_y, overlap_z, overlap, &
!$OMP alpha, beta,i,j,c,d_a_2,d_2,deriv_tmp, &
!$OMP overlap_x0,overlap_y0,overlap_z0) &
!$OMP SHARED(nucl_coord,ao_power,ao_prim_num, &
!$OMP ao_deriv2_x,ao_deriv2_y,ao_deriv2_z,ao_num,ao_coef_normalized_ordered_transp,ao_nucl, &
!$OMP ao_expo_ordered_transp,dim1)
do j=1,ao_num
A_center(1) = nucl_coord( ao_nucl(j), 1 )
A_center(2) = nucl_coord( ao_nucl(j), 2 )
A_center(3) = nucl_coord( ao_nucl(j), 3 )
power_A(1) = ao_power( j, 1 )
power_A(2) = ao_power( j, 2 )
power_A(3) = ao_power( j, 3 )
do i= 1,ao_num
ao_deriv2_x(i,j)= 0.d0
ao_deriv2_y(i,j)= 0.d0
ao_deriv2_z(i,j)= 0.d0
B_center(1) = nucl_coord( ao_nucl(i), 1 )
B_center(2) = nucl_coord( ao_nucl(i), 2 )
B_center(3) = nucl_coord( ao_nucl(i), 3 )
power_B(1) = ao_power( i, 1 )
power_B(2) = ao_power( i, 2 )
power_B(3) = ao_power( i, 3 )
do n = 1,ao_prim_num(j)
alpha = ao_expo_ordered_transp(n,j)
do l = 1, ao_prim_num(i)
beta = ao_expo_ordered_transp(l,i)
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_x0,overlap_y0,overlap_z0,overlap,dim1)
c = ao_coef_normalized_ordered_transp(n,j) * ao_coef_normalized_ordered_transp(l,i)
do j = 1, ao_num
do i = 1, ao_num
ao_deriv2_x(i,j) = ao_deriv2_cosgtos_x(i,j)
ao_deriv2_y(i,j) = ao_deriv2_cosgtos_y(i,j)
ao_deriv2_z(i,j) = ao_deriv2_cosgtos_z(i,j)
enddo
enddo
power_A(1) = power_A(1)-2
if (power_A(1)>-1) then
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,d_a_2,overlap_y,overlap_z,overlap,dim1)
else
d_a_2 = 0.d0
endif
power_A(1) = power_A(1)+4
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,d_2,overlap_y,overlap_z,overlap,dim1)
power_A(1) = power_A(1)-2
else
double precision :: deriv_tmp
deriv_tmp = (-2.d0 * alpha * (2.d0 * power_A(1) +1.d0) * overlap_x0 &
+power_A(1) * (power_A(1)-1.d0) * d_a_2 &
+4.d0 * alpha * alpha * d_2 )*overlap_y0*overlap_z0
dim1=100
ao_deriv2_x(i,j) += c*deriv_tmp
power_A(2) = power_A(2)-2
if (power_A(2)>-1) then
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,d_a_2,overlap_z,overlap,dim1)
else
d_a_2 = 0.d0
endif
power_A(2) = power_A(2)+4
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,d_2,overlap_z,overlap,dim1)
power_A(2) = power_A(2)-2
! -- Dummy call to provide everything
A_center(:) = 0.d0
B_center(:) = 1.d0
alpha = 1.d0
beta = .1d0
power_A = 1
power_B = 0
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,d_a_2,overlap_z,overlap,dim1)
! --
deriv_tmp = (-2.d0 * alpha * (2.d0 * power_A(2) +1.d0 ) * overlap_y0 &
+power_A(2) * (power_A(2)-1.d0) * d_a_2 &
+4.d0 * alpha * alpha * d_2 )*overlap_x0*overlap_z0
ao_deriv2_y(i,j) += c*deriv_tmp
!$OMP PARALLEL DO SCHEDULE(GUIDED) &
!$OMP DEFAULT(NONE) &
!$OMP PRIVATE(A_center,B_center,power_A,power_B,&
!$OMP overlap_y, overlap_z, overlap, &
!$OMP alpha, beta,i,j,c,d_a_2,d_2,deriv_tmp, &
!$OMP overlap_x0,overlap_y0,overlap_z0) &
!$OMP SHARED(nucl_coord,ao_power,ao_prim_num, &
!$OMP ao_deriv2_x,ao_deriv2_y,ao_deriv2_z,ao_num,ao_coef_normalized_ordered_transp,ao_nucl, &
!$OMP ao_expo_ordered_transp,dim1)
do j=1,ao_num
A_center(1) = nucl_coord( ao_nucl(j), 1 )
A_center(2) = nucl_coord( ao_nucl(j), 2 )
A_center(3) = nucl_coord( ao_nucl(j), 3 )
power_A(1) = ao_power( j, 1 )
power_A(2) = ao_power( j, 2 )
power_A(3) = ao_power( j, 3 )
do i= 1,ao_num
ao_deriv2_x(i,j)= 0.d0
ao_deriv2_y(i,j)= 0.d0
ao_deriv2_z(i,j)= 0.d0
B_center(1) = nucl_coord( ao_nucl(i), 1 )
B_center(2) = nucl_coord( ao_nucl(i), 2 )
B_center(3) = nucl_coord( ao_nucl(i), 3 )
power_B(1) = ao_power( i, 1 )
power_B(2) = ao_power( i, 2 )
power_B(3) = ao_power( i, 3 )
do n = 1,ao_prim_num(j)
alpha = ao_expo_ordered_transp(n,j)
do l = 1, ao_prim_num(i)
beta = ao_expo_ordered_transp(l,i)
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_x0,overlap_y0,overlap_z0,overlap,dim1)
c = ao_coef_normalized_ordered_transp(n,j) * ao_coef_normalized_ordered_transp(l,i)
power_A(3) = power_A(3)-2
if (power_A(3)>-1) then
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,overlap_z,d_a_2,overlap,dim1)
else
d_a_2 = 0.d0
endif
power_A(3) = power_A(3)+4
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,overlap_z,d_2,overlap,dim1)
power_A(3) = power_A(3)-2
power_A(1) = power_A(1)-2
if (power_A(1)>-1) then
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,d_a_2,overlap_y,overlap_z,overlap,dim1)
else
d_a_2 = 0.d0
endif
power_A(1) = power_A(1)+4
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,d_2,overlap_y,overlap_z,overlap,dim1)
power_A(1) = power_A(1)-2
deriv_tmp = (-2.d0 * alpha * (2.d0 * power_A(3) +1.d0 ) * overlap_z0 &
+power_A(3) * (power_A(3)-1.d0) * d_a_2 &
+4.d0 * alpha * alpha * d_2 )*overlap_x0*overlap_y0
ao_deriv2_z(i,j) += c*deriv_tmp
double precision :: deriv_tmp
deriv_tmp = (-2.d0 * alpha * (2.d0 * power_A(1) +1.d0) * overlap_x0 &
+power_A(1) * (power_A(1)-1.d0) * d_a_2 &
+4.d0 * alpha * alpha * d_2 )*overlap_y0*overlap_z0
ao_deriv2_x(i,j) += c*deriv_tmp
power_A(2) = power_A(2)-2
if (power_A(2)>-1) then
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,d_a_2,overlap_z,overlap,dim1)
else
d_a_2 = 0.d0
endif
power_A(2) = power_A(2)+4
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,d_2,overlap_z,overlap,dim1)
power_A(2) = power_A(2)-2
deriv_tmp = (-2.d0 * alpha * (2.d0 * power_A(2) +1.d0 ) * overlap_y0 &
+power_A(2) * (power_A(2)-1.d0) * d_a_2 &
+4.d0 * alpha * alpha * d_2 )*overlap_x0*overlap_z0
ao_deriv2_y(i,j) += c*deriv_tmp
power_A(3) = power_A(3)-2
if (power_A(3)>-1) then
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,overlap_z,d_a_2,overlap,dim1)
else
d_a_2 = 0.d0
endif
power_A(3) = power_A(3)+4
call overlap_gaussian_xyz(A_center,B_center,alpha,beta,power_A,power_B,overlap_y,overlap_z,d_2,overlap,dim1)
power_A(3) = power_A(3)-2
deriv_tmp = (-2.d0 * alpha * (2.d0 * power_A(3) +1.d0 ) * overlap_z0 &
+power_A(3) * (power_A(3)-1.d0) * d_a_2 &
+4.d0 * alpha * alpha * d_2 )*overlap_x0*overlap_y0
ao_deriv2_z(i,j) += c*deriv_tmp
enddo
enddo
enddo
enddo
enddo
enddo
!$OMP END PARALLEL DO
!$OMP END PARALLEL DO
endif
END_PROVIDER
! ---
BEGIN_PROVIDER [double precision, ao_kinetic_integrals, (ao_num,ao_num)]
implicit none
BEGIN_DOC

View File

@ -1,3 +1,6 @@
! ---
subroutine give_all_erf_kl_ao(integrals_ao,mu_in,C_center)
implicit none
BEGIN_DOC
@ -15,36 +18,104 @@ subroutine give_all_erf_kl_ao(integrals_ao,mu_in,C_center)
enddo
end
! ---
double precision function NAI_pol_mult_erf_ao(i_ao, j_ao, mu_in, C_center)
double precision function NAI_pol_mult_erf_ao(i_ao,j_ao,mu_in,C_center)
implicit none
BEGIN_DOC
!
! Computes the following integral :
! $\int_{-\infty}^{infty} dr \chi_i(r) \chi_j(r) \frac{\erf(\mu | r - R_C | )}{ | r - R_C | }$.
! $\int_{-\infty}^{infty} dr \chi_i(r) \chi_j(r) \frac{\erf(\mu |r - R_C|)}{|r - R_C|}$.
!
END_DOC
integer, intent(in) :: i_ao,j_ao
implicit none
integer, intent(in) :: i_ao, j_ao
double precision, intent(in) :: mu_in, C_center(3)
integer :: i,j,num_A,num_B, power_A(3), power_B(3), n_pt_in
double precision :: A_center(3), B_center(3),integral, alpha,beta
integer :: i, j, num_A, num_B, power_A(3), power_B(3), n_pt_in
double precision :: A_center(3), B_center(3), integral, alpha, beta
double precision :: NAI_pol_mult_erf
num_A = ao_nucl(i_ao)
power_A(1:3)= ao_power(i_ao,1:3)
num_A = ao_nucl(i_ao)
power_A(1:3) = ao_power(i_ao,1:3)
A_center(1:3) = nucl_coord(num_A,1:3)
num_B = ao_nucl(j_ao)
power_B(1:3)= ao_power(j_ao,1:3)
num_B = ao_nucl(j_ao)
power_B(1:3) = ao_power(j_ao,1:3)
B_center(1:3) = nucl_coord(num_B,1:3)
n_pt_in = n_pt_max_integrals
NAI_pol_mult_erf_ao = 0.d0
do i = 1, ao_prim_num(i_ao)
alpha = ao_expo_ordered_transp(i,i_ao)
do j = 1, ao_prim_num(j_ao)
beta = ao_expo_ordered_transp(j,j_ao)
integral = NAI_pol_mult_erf(A_center,B_center,power_A,power_B,alpha,beta,C_center,n_pt_in,mu_in)
NAI_pol_mult_erf_ao += integral * ao_coef_normalized_ordered_transp(j,j_ao)*ao_coef_normalized_ordered_transp(i,i_ao)
integral = NAI_pol_mult_erf(A_center, B_center, power_A, power_B, alpha, beta, C_center, n_pt_in,mu_in)
NAI_pol_mult_erf_ao += integral * ao_coef_normalized_ordered_transp(j,j_ao) * ao_coef_normalized_ordered_transp(i,i_ao)
enddo
enddo
end
end function NAI_pol_mult_erf_ao
! ---
double precision function NAI_pol_mult_erf_ao_with1s(i_ao, j_ao, beta, B_center, mu_in, C_center)
BEGIN_DOC
!
! Computes the following integral :
! $\int_{-\infty}^{infty} dr \chi_i(r) \chi_j(r) e^{-\beta (r - B_center)^2} \frac{\erf(\mu |r - R_C|)}{|r - R_C|}$.
!
END_DOC
implicit none
integer, intent(in) :: i_ao, j_ao
double precision, intent(in) :: beta, B_center(3)
double precision, intent(in) :: mu_in, C_center(3)
integer :: i, j, power_A1(3), power_A2(3), n_pt_in
double precision :: A1_center(3), A2_center(3), alpha1, alpha2, coef12, coef1, integral
double precision, external :: NAI_pol_mult_erf_with1s, NAI_pol_mult_erf_ao
ASSERT(beta .ge. 0.d0)
if(beta .lt. 1d-10) then
NAI_pol_mult_erf_ao_with1s = NAI_pol_mult_erf_ao(i_ao, j_ao, mu_in, C_center)
return
endif
power_A1(1:3) = ao_power(i_ao,1:3)
power_A2(1:3) = ao_power(j_ao,1:3)
A1_center(1:3) = nucl_coord(ao_nucl(i_ao),1:3)
A2_center(1:3) = nucl_coord(ao_nucl(j_ao),1:3)
n_pt_in = n_pt_max_integrals
NAI_pol_mult_erf_ao_with1s = 0.d0
do i = 1, ao_prim_num(i_ao)
alpha1 = ao_expo_ordered_transp (i,i_ao)
coef1 = ao_coef_normalized_ordered_transp(i,i_ao)
do j = 1, ao_prim_num(j_ao)
alpha2 = ao_expo_ordered_transp(j,j_ao)
coef12 = coef1 * ao_coef_normalized_ordered_transp(j,j_ao)
if(dabs(coef12) .lt. 1d-14) cycle
integral = NAI_pol_mult_erf_with1s( A1_center, A2_center, power_A1, power_A2, alpha1, alpha2 &
, beta, B_center, C_center, n_pt_in, mu_in )
NAI_pol_mult_erf_ao_with1s += integral * coef12
enddo
enddo
end function NAI_pol_mult_erf_ao_with1s
! ---
double precision function NAI_pol_mult_erf(A_center, B_center, power_A, power_B, alpha, beta, C_center, n_pt_in, mu_in)
@ -127,58 +198,221 @@ end function NAI_pol_mult_erf
! ---
double precision function NAI_pol_mult_erf_ao_with1s(i_ao, j_ao, beta, B_center, mu_in, C_center)
subroutine NAI_pol_mult_erf_v(A_center, B_center, power_A, power_B, alpha, beta, C_center, LD_C, n_pt_in, mu_in, res_v, LD_resv, n_points)
BEGIN_DOC
!
! Computes the following integral :
! $\int_{-\infty}^{infty} dr \chi_i(r) \chi_j(r) e^{-\beta (r - B_center)^2} \frac{\erf(\mu |r - R_C|)}{|r - R_C|}$.
!
! .. math::
!
! \int dr (x-A_x)^a (x-B_x)^b \exp(-\alpha (x-A_x)^2 - \beta (x-B_x)^2 )
! \frac{\erf(\mu |r - R_C |)}{| r - R_C |}$.
!
END_DOC
include 'utils/constants.include.F'
implicit none
integer, intent(in) :: i_ao, j_ao
double precision, intent(in) :: beta, B_center(3)
double precision, intent(in) :: mu_in, C_center(3)
integer :: i, j, power_A1(3), power_A2(3), n_pt_in
double precision :: A1_center(3), A2_center(3), alpha1, alpha2, coef12, coef1, integral
integer, intent(in) :: n_pt_in, n_points, LD_C, LD_resv
integer, intent(in) :: power_A(3), power_B(3)
double precision, intent(in) :: A_center(3), B_center(3), alpha, beta, mu_in
double precision, intent(in) :: C_center(LD_C,3)
double precision, intent(out) :: res_v(LD_resv)
double precision, external :: NAI_pol_mult_erf_with1s, NAI_pol_mult_erf_ao
integer :: i, n_pt, n_pt_out, ipoint
double precision :: P_center(3)
double precision :: d(0:n_pt_in), coeff, dist, const, factor
double precision :: const_factor, dist_integral
double precision :: accu, p_inv, p, rho, p_inv_2
double precision :: p_new, p_new2, coef_tmp
ASSERT(beta .ge. 0.d0)
if(beta .lt. 1d-10) then
NAI_pol_mult_erf_ao_with1s = NAI_pol_mult_erf_ao(i_ao, j_ao, mu_in, C_center)
double precision :: rint
res_V(1:LD_resv) = 0.d0
p = alpha + beta
p_inv = 1.d0 / p
p_inv_2 = 0.5d0 * p_inv
rho = alpha * beta * p_inv
p_new = mu_in / dsqrt(p + mu_in * mu_in)
p_new2 = p_new * p_new
coef_tmp = p * p_new2
dist = 0.d0
do i = 1, 3
P_center(i) = (alpha * A_center(i) + beta * B_center(i)) * p_inv
dist += (A_center(i) - B_center(i)) * (A_center(i) - B_center(i))
enddo
const_factor = dist * rho
if(const_factor > 80.d0) then
return
endif
factor = dexp(-const_factor)
coeff = dtwo_pi * factor * p_inv * p_new
n_pt = 2 * ( power_A(1) + power_B(1) + power_A(2) + power_B(2) + power_A(3) + power_B(3) )
if(n_pt == 0) then
do ipoint = 1, n_points
dist_integral = 0.d0
do i = 1, 3
dist_integral += (P_center(i) - C_center(ipoint,i)) * (P_center(i) - C_center(ipoint,i))
enddo
const = coef_tmp * dist_integral
res_v(ipoint) = coeff * rint(0, const)
enddo
else
do ipoint = 1, n_points
dist_integral = 0.d0
do i = 1, 3
dist_integral += (P_center(i) - C_center(ipoint,i)) * (P_center(i) - C_center(ipoint,i))
enddo
const = coef_tmp * dist_integral
do i = 0, n_pt_in
d(i) = 0.d0
enddo
call give_polynomial_mult_center_one_e_erf_opt(A_center, B_center, power_A, power_B, C_center(ipoint,1:3), n_pt_in, d, n_pt_out, p_inv_2, p_new2, P_center)
if(n_pt_out < 0) then
res_v(ipoint) = 0.d0
cycle
endif
! sum of integrals of type : int {t,[0,1]} exp-(rho.(P-Q)^2 * t^2) * t^i
accu = 0.d0
do i = 0, n_pt_out, 2
accu += d(i) * rint(i/2, const)
enddo
res_v(ipoint) = accu * coeff
enddo
endif
end subroutine NAI_pol_mult_erf_v
! ---
double precision function NAI_pol_mult_erf_with1s( A1_center, A2_center, power_A1, power_A2, alpha1, alpha2 &
, beta, B_center, C_center, n_pt_in, mu_in )
BEGIN_DOC
!
! Computes the following integral :
!
! .. math::
!
! \int dx (x - A1_x)^a_1 (x - B1_x)^a_2 \exp(-\alpha_1 (x - A1_x)^2 - \alpha_2 (x - A2_x)^2)
! \int dy (y - A1_y)^b_1 (y - B1_y)^b_2 \exp(-\alpha_1 (y - A1_y)^2 - \alpha_2 (y - A2_y)^2)
! \int dz (x - A1_z)^c_1 (z - B1_z)^c_2 \exp(-\alpha_1 (z - A1_z)^2 - \alpha_2 (z - A2_z)^2)
! \exp(-\beta (r - B)^2)
! \frac{\erf(\mu |r - R_C|)}{|r - R_C|}$.
!
END_DOC
include 'utils/constants.include.F'
implicit none
integer, intent(in) :: n_pt_in
integer, intent(in) :: power_A1(3), power_A2(3)
double precision, intent(in) :: C_center(3), A1_center(3), A2_center(3), B_center(3)
double precision, intent(in) :: alpha1, alpha2, beta, mu_in
integer :: i, n_pt, n_pt_out
double precision :: alpha12, alpha12_inv, alpha12_inv_2, rho12, A12_center(3), dist12, const_factor12
double precision :: p, p_inv, p_inv_2, rho, P_center(3), dist, const_factor
double precision :: dist_integral
double precision :: d(0:n_pt_in), coeff, const, factor
double precision :: accu
double precision :: p_new
double precision :: rint
! e^{-alpha1 (r - A1)^2} e^{-alpha2 (r - A2)^2} = e^{-K12} e^{-alpha12 (r - A12)^2}
alpha12 = alpha1 + alpha2
alpha12_inv = 1.d0 / alpha12
alpha12_inv_2 = 0.5d0 * alpha12_inv
rho12 = alpha1 * alpha2 * alpha12_inv
A12_center(1) = (alpha1 * A1_center(1) + alpha2 * A2_center(1)) * alpha12_inv
A12_center(2) = (alpha1 * A1_center(2) + alpha2 * A2_center(2)) * alpha12_inv
A12_center(3) = (alpha1 * A1_center(3) + alpha2 * A2_center(3)) * alpha12_inv
dist12 = (A1_center(1) - A2_center(1)) * (A1_center(1) - A2_center(1)) &
+ (A1_center(2) - A2_center(2)) * (A1_center(2) - A2_center(2)) &
+ (A1_center(3) - A2_center(3)) * (A1_center(3) - A2_center(3))
const_factor12 = dist12 * rho12
if(const_factor12 > 80.d0) then
NAI_pol_mult_erf_with1s = 0.d0
return
endif
power_A1(1:3) = ao_power(i_ao,1:3)
power_A2(1:3) = ao_power(j_ao,1:3)
! ---
A1_center(1:3) = nucl_coord(ao_nucl(i_ao),1:3)
A2_center(1:3) = nucl_coord(ao_nucl(j_ao),1:3)
! e^{-K12} e^{-alpha12 (r - A12)^2} e^{-beta (r - B)^2} = e^{-K} e^{-p (r - P)^2}
p = alpha12 + beta
p_inv = 1.d0 / p
p_inv_2 = 0.5d0 * p_inv
rho = alpha12 * beta * p_inv
P_center(1) = (alpha12 * A12_center(1) + beta * B_center(1)) * p_inv
P_center(2) = (alpha12 * A12_center(2) + beta * B_center(2)) * p_inv
P_center(3) = (alpha12 * A12_center(3) + beta * B_center(3)) * p_inv
dist = (A12_center(1) - B_center(1)) * (A12_center(1) - B_center(1)) &
+ (A12_center(2) - B_center(2)) * (A12_center(2) - B_center(2)) &
+ (A12_center(3) - B_center(3)) * (A12_center(3) - B_center(3))
n_pt_in = n_pt_max_integrals
const_factor = const_factor12 + dist * rho
if(const_factor > 80.d0) then
NAI_pol_mult_erf_with1s = 0.d0
return
endif
NAI_pol_mult_erf_ao_with1s = 0.d0
do i = 1, ao_prim_num(i_ao)
alpha1 = ao_expo_ordered_transp (i,i_ao)
coef1 = ao_coef_normalized_ordered_transp(i,i_ao)
dist_integral = (P_center(1) - C_center(1)) * (P_center(1) - C_center(1)) &
+ (P_center(2) - C_center(2)) * (P_center(2) - C_center(2)) &
+ (P_center(3) - C_center(3)) * (P_center(3) - C_center(3))
do j = 1, ao_prim_num(j_ao)
alpha2 = ao_expo_ordered_transp(j,j_ao)
coef12 = coef1 * ao_coef_normalized_ordered_transp(j,j_ao)
if(dabs(coef12) .lt. 1d-14) cycle
! ---
integral = NAI_pol_mult_erf_with1s( A1_center, A2_center, power_A1, power_A2, alpha1, alpha2 &
, beta, B_center, C_center, n_pt_in, mu_in )
p_new = mu_in / dsqrt(p + mu_in * mu_in)
factor = dexp(-const_factor)
coeff = dtwo_pi * factor * p_inv * p_new
NAI_pol_mult_erf_ao_with1s += integral * coef12
enddo
n_pt = 2 * ( (power_A1(1) + power_A2(1)) + (power_A1(2) + power_A2(2)) + (power_A1(3) + power_A2(3)) )
const = p * dist_integral * p_new * p_new
if(n_pt == 0) then
NAI_pol_mult_erf_with1s = coeff * rint(0, const)
return
endif
do i = 0, n_pt_in
d(i) = 0.d0
enddo
p_new = p_new * p_new
call give_polynomial_mult_center_one_e_erf_opt( A1_center, A2_center, power_A1, power_A2, C_center, n_pt_in, d, n_pt_out, p_inv_2, p_new, P_center)
end function NAI_pol_mult_erf_ao_with1s
if(n_pt_out < 0) then
NAI_pol_mult_erf_with1s = 0.d0
return
endif
! sum of integrals of type : int {t,[0,1]} exp-(rho.(P-Q)^2 * t^2) * t^i
accu = 0.d0
do i = 0, n_pt_out, 2
accu += d(i) * rint(i/2, const)
enddo
NAI_pol_mult_erf_with1s = accu * coeff
end function NAI_pol_mult_erf_with1s
! ---
subroutine NAI_pol_mult_erf_with1s_v(A1_center, A2_center, power_A1, power_A2, alpha1, alpha2, beta, B_center, LD_B, C_center, LD_C, n_pt_in, mu_in, res_v, LD_resv, n_points)
@ -428,107 +662,6 @@ subroutine give_polynomial_mult_center_one_e_erf_opt(A_center, B_center, power_A
end subroutine give_polynomial_mult_center_one_e_erf_opt
! ---
subroutine NAI_pol_mult_erf_v(A_center, B_center, power_A, power_B, alpha, beta, C_center, LD_C, n_pt_in, mu_in, res_v, LD_resv, n_points)
BEGIN_DOC
!
! Computes the following integral :
!
! .. math::
!
! \int dr (x-A_x)^a (x-B_x)^b \exp(-\alpha (x-A_x)^2 - \beta (x-B_x)^2 )
! \frac{\erf(\mu |r - R_C |)}{| r - R_C |}$.
!
END_DOC
include 'utils/constants.include.F'
implicit none
integer, intent(in) :: n_pt_in, n_points, LD_C, LD_resv
integer, intent(in) :: power_A(3), power_B(3)
double precision, intent(in) :: A_center(3), B_center(3), alpha, beta, mu_in
double precision, intent(in) :: C_center(LD_C,3)
double precision, intent(out) :: res_v(LD_resv)
integer :: i, n_pt, n_pt_out, ipoint
double precision :: P_center(3)
double precision :: d(0:n_pt_in), coeff, dist, const, factor
double precision :: const_factor, dist_integral
double precision :: accu, p_inv, p, rho, p_inv_2
double precision :: p_new, p_new2, coef_tmp
double precision :: rint
res_V(1:LD_resv) = 0.d0
p = alpha + beta
p_inv = 1.d0 / p
p_inv_2 = 0.5d0 * p_inv
rho = alpha * beta * p_inv
p_new = mu_in / dsqrt(p + mu_in * mu_in)
p_new2 = p_new * p_new
coef_tmp = p * p_new2
dist = 0.d0
do i = 1, 3
P_center(i) = (alpha * A_center(i) + beta * B_center(i)) * p_inv
dist += (A_center(i) - B_center(i)) * (A_center(i) - B_center(i))
enddo
const_factor = dist * rho
if(const_factor > 80.d0) then
return
endif
factor = dexp(-const_factor)
coeff = dtwo_pi * factor * p_inv * p_new
n_pt = 2 * ( power_A(1) + power_B(1) + power_A(2) + power_B(2) + power_A(3) + power_B(3) )
if(n_pt == 0) then
do ipoint = 1, n_points
dist_integral = 0.d0
do i = 1, 3
dist_integral += (P_center(i) - C_center(ipoint,i)) * (P_center(i) - C_center(ipoint,i))
enddo
const = coef_tmp * dist_integral
res_v(ipoint) = coeff * rint(0, const)
enddo
else
do ipoint = 1, n_points
dist_integral = 0.d0
do i = 1, 3
dist_integral += (P_center(i) - C_center(ipoint,i)) * (P_center(i) - C_center(ipoint,i))
enddo
const = coef_tmp * dist_integral
do i = 0, n_pt_in
d(i) = 0.d0
enddo
call give_polynomial_mult_center_one_e_erf_opt(A_center, B_center, power_A, power_B, C_center(ipoint,1:3), n_pt_in, d, n_pt_out, p_inv_2, p_new2, P_center)
if(n_pt_out < 0) then
res_v(ipoint) = 0.d0
cycle
endif
! sum of integrals of type : int {t,[0,1]} exp-(rho.(P-Q)^2 * t^2) * t^i
accu = 0.d0
do i = 0, n_pt_out, 2
accu += d(i) * rint(i/2, const)
enddo
res_v(ipoint) = accu * coeff
enddo
endif
end subroutine NAI_pol_mult_erf_v
subroutine give_polynomial_mult_center_one_e_erf(A_center,B_center,alpha,beta,power_A,power_B,C_center,n_pt_in,d,n_pt_out,mu_in)
@ -659,113 +792,3 @@ subroutine give_polynomial_mult_center_one_e_erf(A_center,B_center,alpha,beta,po
end
double precision function NAI_pol_mult_erf_with1s( A1_center, A2_center, power_A1, power_A2, alpha1, alpha2 &
, beta, B_center, C_center, n_pt_in, mu_in )
BEGIN_DOC
!
! Computes the following integral :
!
! .. math::
!
! \int dx (x - A1_x)^a_1 (x - B1_x)^a_2 \exp(-\alpha_1 (x - A1_x)^2 - \alpha_2 (x - A2_x)^2)
! \int dy (y - A1_y)^b_1 (y - B1_y)^b_2 \exp(-\alpha_1 (y - A1_y)^2 - \alpha_2 (y - A2_y)^2)
! \int dz (x - A1_z)^c_1 (z - B1_z)^c_2 \exp(-\alpha_1 (z - A1_z)^2 - \alpha_2 (z - A2_z)^2)
! \exp(-\beta (r - B)^2)
! \frac{\erf(\mu |r - R_C|)}{|r - R_C|}$.
!
END_DOC
include 'utils/constants.include.F'
implicit none
integer, intent(in) :: n_pt_in
integer, intent(in) :: power_A1(3), power_A2(3)
double precision, intent(in) :: C_center(3), A1_center(3), A2_center(3), B_center(3)
double precision, intent(in) :: alpha1, alpha2, beta, mu_in
integer :: i, n_pt, n_pt_out
double precision :: alpha12, alpha12_inv, alpha12_inv_2, rho12, A12_center(3), dist12, const_factor12
double precision :: p, p_inv, p_inv_2, rho, P_center(3), dist, const_factor
double precision :: dist_integral
double precision :: d(0:n_pt_in), coeff, const, factor
double precision :: accu
double precision :: p_new
double precision :: rint
! e^{-alpha1 (r - A1)^2} e^{-alpha2 (r - A2)^2} = e^{-K12} e^{-alpha12 (r - A12)^2}
alpha12 = alpha1 + alpha2
alpha12_inv = 1.d0 / alpha12
alpha12_inv_2 = 0.5d0 * alpha12_inv
rho12 = alpha1 * alpha2 * alpha12_inv
A12_center(1) = (alpha1 * A1_center(1) + alpha2 * A2_center(1)) * alpha12_inv
A12_center(2) = (alpha1 * A1_center(2) + alpha2 * A2_center(2)) * alpha12_inv
A12_center(3) = (alpha1 * A1_center(3) + alpha2 * A2_center(3)) * alpha12_inv
dist12 = (A1_center(1) - A2_center(1)) * (A1_center(1) - A2_center(1)) &
+ (A1_center(2) - A2_center(2)) * (A1_center(2) - A2_center(2)) &
+ (A1_center(3) - A2_center(3)) * (A1_center(3) - A2_center(3))
const_factor12 = dist12 * rho12
if(const_factor12 > 80.d0) then
NAI_pol_mult_erf_with1s = 0.d0
return
endif
! ---
! e^{-K12} e^{-alpha12 (r - A12)^2} e^{-beta (r - B)^2} = e^{-K} e^{-p (r - P)^2}
p = alpha12 + beta
p_inv = 1.d0 / p
p_inv_2 = 0.5d0 * p_inv
rho = alpha12 * beta * p_inv
P_center(1) = (alpha12 * A12_center(1) + beta * B_center(1)) * p_inv
P_center(2) = (alpha12 * A12_center(2) + beta * B_center(2)) * p_inv
P_center(3) = (alpha12 * A12_center(3) + beta * B_center(3)) * p_inv
dist = (A12_center(1) - B_center(1)) * (A12_center(1) - B_center(1)) &
+ (A12_center(2) - B_center(2)) * (A12_center(2) - B_center(2)) &
+ (A12_center(3) - B_center(3)) * (A12_center(3) - B_center(3))
const_factor = const_factor12 + dist * rho
if(const_factor > 80.d0) then
NAI_pol_mult_erf_with1s = 0.d0
return
endif
dist_integral = (P_center(1) - C_center(1)) * (P_center(1) - C_center(1)) &
+ (P_center(2) - C_center(2)) * (P_center(2) - C_center(2)) &
+ (P_center(3) - C_center(3)) * (P_center(3) - C_center(3))
! ---
p_new = mu_in / dsqrt(p + mu_in * mu_in)
factor = dexp(-const_factor)
coeff = dtwo_pi * factor * p_inv * p_new
n_pt = 2 * ( (power_A1(1) + power_A2(1)) + (power_A1(2) + power_A2(2)) + (power_A1(3) + power_A2(3)) )
const = p * dist_integral * p_new * p_new
if(n_pt == 0) then
NAI_pol_mult_erf_with1s = coeff * rint(0, const)
return
endif
do i = 0, n_pt_in
d(i) = 0.d0
enddo
p_new = p_new * p_new
call give_polynomial_mult_center_one_e_erf_opt( A1_center, A2_center, power_A1, power_A2, C_center, n_pt_in, d, n_pt_out, p_inv_2, p_new, P_center)
if(n_pt_out < 0) then
NAI_pol_mult_erf_with1s = 0.d0
return
endif
! sum of integrals of type : int {t,[0,1]} exp-(rho.(P-Q)^2 * t^2) * t^i
accu = 0.d0
do i = 0, n_pt_out, 2
accu += d(i) * rint(i/2, const)
enddo
NAI_pol_mult_erf_with1s = accu * coeff
end function NAI_pol_mult_erf_with1s

View File

<
@ -1,4 +1,8 @@
! ---
BEGIN_PROVIDER [ double precision, ao_integrals_n_e, (ao_num,ao_num)]
BEGIN_DOC
! Nucleus-electron interaction, in the |AO| basis set.
!
@ -6,84 +10,100 @@ BEGIN_PROVIDER [ double precision, ao_integrals_n_e, (ao_num,ao_num)]
!
! These integrals also contain the pseudopotential integrals.
END_DOC
implicit none
double precision :: alpha, beta, gama, delta
integer :: num_A,num_B
double precision :: A_center(3),B_center(3),C_center(3)
integer :: power_A(3),power_B(3)
integer :: i,j,k,l,n_pt_in,m
double precision :: overlap_x,overlap_y,overlap_z,overlap,dx,NAI_pol_mult
integer :: num_A, num_B, power_A(3), power_B(3)
integer :: i, j, k, l, n_pt_in, m
double precision :: alpha, beta
double precision :: A_center(3),B_center(3),C_center(3)
double precision :: overlap_x,overlap_y,overlap_z,overlap,dx,NAI_pol_mult
ao_integrals_n_e = 0.d0
if (read_ao_integrals_n_e) then
call ezfio_get_ao_one_e_ints_ao_integrals_n_e(ao_integrals_n_e)
print *, 'AO N-e integrals read from disk'
else
ao_integrals_n_e = 0.d0
if(use_cosgtos) then
!print *, " use_cosgtos for ao_integrals_n_e ?", use_cosgtos
! _
! /| / |_)
! | / | \
!
do j = 1, ao_num
do i = 1, ao_num
ao_integrals_n_e(i,j) = ao_integrals_n_e_cosgtos(i,j)
enddo
enddo
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (i,j,k,l,m,alpha,beta,A_center,B_center,C_center,power_A,power_B,&
!$OMP num_A,num_B,Z,c,n_pt_in) &
!$OMP SHARED (ao_num,ao_prim_num,ao_expo_ordered_transp,ao_power,ao_nucl,nucl_coord,ao_coef_normalized_ordered_transp,&
!$OMP n_pt_max_integrals,ao_integrals_n_e,nucl_num,nucl_charge)
else