9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-22 19:43:32 +01:00

Change 1.d-32 into 1.d-20 in pseudo

This commit is contained in:
Anthony Scemama 2022-06-10 11:34:40 +02:00
parent 5c66b316b9
commit a012a078cd
4 changed files with 21 additions and 19 deletions

@ -1 +1 @@
Subproject commit 242151e03d1d6bf042387226431d82d35845686a Subproject commit 90ee61f5041c7c94a0c605625a264860292813a0

View File

@ -1869,7 +1869,7 @@ double precision function int_prod_bessel(l,gam,n,m,a,b,arg)
qk = dble(q) qk = dble(q)
two_qkmp1 = 2.d0*(qk+mk)+1.d0 two_qkmp1 = 2.d0*(qk+mk)+1.d0
do k=0,q-1 do k=0,q-1
if (s_q_k < 1.d-32) then if (s_q_k < 1.d-20) then
s_q_k = 0.d0 s_q_k = 0.d0
exit exit
endif endif

View File

@ -142,7 +142,7 @@ subroutine ao_idx2_sq(i,j,ij)
ij=i*i ij=i*i
endif endif
end end
subroutine idx2_tri_int(i,j,ij) subroutine idx2_tri_int(i,j,ij)
implicit none implicit none
integer, intent(in) :: i,j integer, intent(in) :: i,j
@ -152,7 +152,7 @@ subroutine idx2_tri_int(i,j,ij)
q = min(i,j) q = min(i,j)
ij = q+ishft(p*p-p,-1) ij = q+ishft(p*p-p,-1)
end end
subroutine ao_idx2_tri_key(i,j,ij) subroutine ao_idx2_tri_key(i,j,ij)
use map_module use map_module
implicit none implicit none
@ -163,8 +163,8 @@ subroutine ao_idx2_tri_key(i,j,ij)
q = min(i,j) q = min(i,j)
ij = q+ishft(p*p-p,-1) ij = q+ishft(p*p-p,-1)
end end
subroutine two_e_integrals_index_2fold(i,j,k,l,i1) subroutine two_e_integrals_index_2fold(i,j,k,l,i1)
use map_module use map_module
implicit none implicit none
integer, intent(in) :: i,j,k,l integer, intent(in) :: i,j,k,l
@ -176,7 +176,7 @@ subroutine two_e_integrals_index_2fold(i,j,k,l,i1)
call ao_idx2_tri_key(ik,jl,i1) call ao_idx2_tri_key(ik,jl,i1)
end end
subroutine ao_idx2_sq_rev(i,k,ik) subroutine ao_idx2_sq_rev(i,k,ik)
BEGIN_DOC BEGIN_DOC
! reverse square compound index ! reverse square compound index
END_DOC END_DOC
@ -398,7 +398,7 @@ BEGIN_PROVIDER [ complex*16, ao_integrals_cache_periodic, (0:64*64*64*64) ]
tmp_im = 0.d0 tmp_im = 0.d0
integral = dcmplx(tmp_re,tmp_im) integral = dcmplx(tmp_re,tmp_im)
endif endif
ii = l-ao_integrals_cache_min ii = l-ao_integrals_cache_min
ii = ior( shiftl(ii,6), k-ao_integrals_cache_min) ii = ior( shiftl(ii,6), k-ao_integrals_cache_min)
ii = ior( shiftl(ii,6), j-ao_integrals_cache_min) ii = ior( shiftl(ii,6), j-ao_integrals_cache_min)
@ -473,7 +473,7 @@ subroutine get_ao_two_e_integrals(j,k,l,sze,out_val)
BEGIN_DOC BEGIN_DOC
! Gets multiple AO bi-electronic integral from the AO map . ! Gets multiple AO bi-electronic integral from the AO map .
! All i are retrieved for j,k,l fixed. ! All i are retrieved for j,k,l fixed.
! physicist convention : <ij|kl> ! physicist convention : <ij|kl>
END_DOC END_DOC
implicit none implicit none
integer, intent(in) :: j,k,l, sze integer, intent(in) :: j,k,l, sze
@ -482,7 +482,7 @@ subroutine get_ao_two_e_integrals(j,k,l,sze,out_val)
integer :: i integer :: i
integer(key_kind) :: hash integer(key_kind) :: hash
logical, external :: ao_one_e_integral_zero logical, external :: ao_one_e_integral_zero
PROVIDE ao_two_e_integrals_in_map ao_integrals_map PROVIDE ao_two_e_integrals_in_map ao_integrals_map
if (ao_one_e_integral_zero(j,l)) then if (ao_one_e_integral_zero(j,l)) then
out_val = 0.d0 out_val = 0.d0
@ -502,7 +502,7 @@ subroutine get_ao_two_e_integrals_periodic(j,k,l,sze,out_val)
BEGIN_DOC BEGIN_DOC
! Gets multiple AO bi-electronic integral from the AO map . ! Gets multiple AO bi-electronic integral from the AO map .
! All i are retrieved for j,k,l fixed. ! All i are retrieved for j,k,l fixed.
! physicist convention : <ij|kl> ! physicist convention : <ij|kl>
END_DOC END_DOC
implicit none implicit none
integer, intent(in) :: j,k,l, sze integer, intent(in) :: j,k,l, sze

View File

@ -38,7 +38,7 @@ BEGIN_PROVIDER [ logical, mo_two_e_integrals_in_map ]
print*, 'MO integrals provided' print*, 'MO integrals provided'
return return
else else
PROVIDE ao_two_e_integrals_in_map PROVIDE ao_two_e_integrals_in_map
endif endif
print *, '' print *, ''
@ -53,7 +53,7 @@ BEGIN_PROVIDER [ logical, mo_two_e_integrals_in_map ]
! call four_idx_novvvv ! call four_idx_novvvv
call four_idx_novvvv_old call four_idx_novvvv_old
else else
if (ao_num*ao_num*ao_num*ao_num*32.d-6 < dble(qp_max_mem)) then if (dble(ao_num)**4 * 32.d-9 < dble(qp_max_mem)) then
call four_idx_dgemm call four_idx_dgemm
else else
call add_integrals_to_map(full_ijkl_bitmask_4) call add_integrals_to_map(full_ijkl_bitmask_4)
@ -247,14 +247,14 @@ subroutine add_integrals_to_map(mask_ijkl)
return return
endif endif
size_buffer = min(ao_num*ao_num*ao_num,16000000)
print*, 'Buffers : ', 8.*(mo_num*(n_j)*(n_k+1) + mo_num+&
ao_num+ao_num*ao_num+ size_buffer*3)/(1024*1024), 'MB / core'
double precision :: accu_bis double precision :: accu_bis
accu_bis = 0.d0 accu_bis = 0.d0
call wall_time(wall_1) call wall_time(wall_1)
size_buffer = min( (qp_max_mem/(nproc*5)),mo_num*mo_num*mo_num)
print*, 'Buffers : ', 8.*(mo_num*(n_j)*(n_k+1) + mo_num+&
ao_num+ao_num*ao_num+ size_buffer*3)/(1024*1024), 'MB / core'
!$OMP PARALLEL PRIVATE(l1,k1,j1,i1,i2,i3,i4,i,j,k,l,c, ii1,kmax, & !$OMP PARALLEL PRIVATE(l1,k1,j1,i1,i2,i3,i4,i,j,k,l,c, ii1,kmax, &
!$OMP two_e_tmp_0_idx, two_e_tmp_0, two_e_tmp_1,two_e_tmp_2,two_e_tmp_3,& !$OMP two_e_tmp_0_idx, two_e_tmp_0, two_e_tmp_1,two_e_tmp_2,two_e_tmp_3,&
!$OMP buffer_i,buffer_value,n_integrals,wall_2,i0,j0,k0,l0, & !$OMP buffer_i,buffer_value,n_integrals,wall_2,i0,j0,k0,l0, &
@ -265,6 +265,10 @@ subroutine add_integrals_to_map(mask_ijkl)
!$OMP mo_coef_transp_is_built, list_ijkl, & !$OMP mo_coef_transp_is_built, list_ijkl, &
!$OMP mo_coef_is_built, wall_1, & !$OMP mo_coef_is_built, wall_1, &
!$OMP mo_coef,mo_integrals_threshold,mo_integrals_map) !$OMP mo_coef,mo_integrals_threshold,mo_integrals_map)
thread_num = 0
!$ thread_num = omp_get_thread_num()
n_integrals = 0 n_integrals = 0
wall_0 = wall_1 wall_0 = wall_1
allocate(two_e_tmp_3(mo_num, n_j, n_k), & allocate(two_e_tmp_3(mo_num, n_j, n_k), &
@ -275,8 +279,6 @@ subroutine add_integrals_to_map(mask_ijkl)
buffer_i(size_buffer), & buffer_i(size_buffer), &
buffer_value(size_buffer) ) buffer_value(size_buffer) )
thread_num = 0
!$ thread_num = omp_get_thread_num()
!$OMP DO SCHEDULE(guided) !$OMP DO SCHEDULE(guided)
do l1 = 1,ao_num do l1 = 1,ao_num
two_e_tmp_3 = 0.d0 two_e_tmp_3 = 0.d0