mirror of
https://github.com/LCPQ/quantum_package
synced 2024-12-22 20:35:19 +01:00
Merge branch 'master' of github.com:scemama/quantum_package
This commit is contained in:
commit
45f75131f6
@ -84,8 +84,8 @@ program fci_zmq
|
||||
|
||||
if(do_pt2_end)then
|
||||
print*,'Last iteration only to compute the PT2'
|
||||
threshold_selectors = 1.d0
|
||||
threshold_generators = 0.9999d0
|
||||
threshold_selectors = threshold_selectors_pt2
|
||||
threshold_generators = threshold_generators_pt2
|
||||
TOUCH threshold_selectors threshold_generators
|
||||
E_CI_before(1:N_states) = CI_energy(1:N_states)
|
||||
call ZMQ_selection(0, pt2)
|
||||
|
@ -643,13 +643,12 @@ subroutine davidson_diag_hjj_sjj_mrcc(dets_in,u_in,H_jj,S2_jj,energies,dim_in,sz
|
||||
integer :: i,j,k,l,m
|
||||
logical :: converged
|
||||
|
||||
double precision, allocatable :: overlap(:,:)
|
||||
double precision :: u_dot_v, u_dot_u
|
||||
|
||||
integer :: k_pairs, kl
|
||||
|
||||
integer :: iter2
|
||||
double precision, allocatable :: W(:,:), U(:,:), R(:,:), S(:,:)
|
||||
double precision, allocatable :: W(:,:), U(:,:), S(:,:)
|
||||
double precision, allocatable :: y(:,:), h(:,:), lambda(:), s2(:)
|
||||
double precision, allocatable :: c(:), s_(:,:), s_tmp(:,:)
|
||||
double precision :: diag_h_mat_elem
|
||||
@ -660,7 +659,7 @@ subroutine davidson_diag_hjj_sjj_mrcc(dets_in,u_in,H_jj,S2_jj,energies,dim_in,sz
|
||||
integer :: shift, shift2, itermax
|
||||
include 'constants.include.F'
|
||||
|
||||
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: U, W, R, S, y, h, lambda
|
||||
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: U, W, S, y, h, lambda
|
||||
if (N_st_diag > sze) then
|
||||
stop 'error in Davidson : N_st_diag > sze'
|
||||
endif
|
||||
@ -905,8 +904,8 @@ subroutine davidson_diag_hjj_sjj_mrcc(dets_in,u_in,H_jj,S2_jj,energies,dim_in,sz
|
||||
|
||||
deallocate ( &
|
||||
W, residual_norm, &
|
||||
U, overlap, &
|
||||
R, c, S, &
|
||||
U, &
|
||||
c, S, &
|
||||
h, &
|
||||
y, s_, s_tmp, &
|
||||
lambda &
|
||||
|
@ -229,13 +229,14 @@ subroutine davidson_diag_hjj_sjj(dets_in,u_in,H_jj,S2_jj,energies,dim_in,sze,N_s
|
||||
! Compute h_kl = <u_k | W_l> = <u_k| H |u_l>
|
||||
! -------------------------------------------
|
||||
|
||||
call dgemm('T','N', shift2, N_st_diag, sze, &
|
||||
1.d0, U, size(U,1), W(1,shift+1), size(W,1), &
|
||||
0.d0, h(1,shift+1), size(h,1))
|
||||
call dgemm('T','N', shift2, shift2, sze, &
|
||||
1.d0, U(1,1), size(U,1), W(1,1), size(W,1), &
|
||||
0.d0, h(1,1), size(h,1))
|
||||
|
||||
call dgemm('T','N', shift2, shift2, sze, &
|
||||
1.d0, U(1,1), size(U,1), S(1,1), size(S,1), &
|
||||
0.d0, s_(1,1), size(s_,1))
|
||||
|
||||
call dgemm('T','N', shift2, N_st_diag, sze, &
|
||||
1.d0, U, size(U,1), S(1,shift+1), size(S,1), &
|
||||
0.d0, s_(1,shift+1), size(s_,1))
|
||||
|
||||
! Diagonalize h
|
||||
! -------------
|
||||
|
@ -35,7 +35,7 @@ subroutine occ_pattern_to_dets_size(o,sze,n_alpha,Nint)
|
||||
bmax += popcnt( o(k,1) )
|
||||
amax -= popcnt( o(k,2) )
|
||||
enddo
|
||||
sze = int( min(binom_func(bmax, amax), 1.d8) )
|
||||
sze = 2*int( min(binom_func(bmax, amax), 1.d8) )
|
||||
|
||||
end
|
||||
|
||||
@ -76,27 +76,6 @@ subroutine occ_pattern_to_dets(o,d,sze,n_alpha,Nint)
|
||||
enddo
|
||||
enddo
|
||||
|
||||
! !TODO DEBUG
|
||||
! integer :: j,s
|
||||
! do i=1,nd
|
||||
! do j=1,i-1
|
||||
! na=0
|
||||
! do k=1,Nint
|
||||
! if((d(k,1,j) /= d(k,1,i)).or. &
|
||||
! (d(k,2,j) /= d(k,2,i))) then
|
||||
! s=1
|
||||
! exit
|
||||
! endif
|
||||
! enddo
|
||||
! if ( j== 0 ) then
|
||||
! print *, 'det ',i,' and ',j,' equal:'
|
||||
! call debug_det(d(1,1,j),Nint)
|
||||
! call debug_det(d(1,1,i),Nint)
|
||||
! stop
|
||||
! endif
|
||||
! enddo
|
||||
! enddo
|
||||
! !TODO DEBUG
|
||||
end
|
||||
|
||||
recursive subroutine rec_occ_pattern_to_dets(list_todo,nt,list_a,na,d,nd,sze,amax,Nint)
|
||||
@ -226,26 +205,7 @@ end
|
||||
enddo
|
||||
|
||||
deallocate(iorder,duplicate,bit_tmp,tmp_array)
|
||||
! !TODO DEBUG
|
||||
! integer :: s
|
||||
! do i=1,N_occ_pattern
|
||||
! do j=i+1,N_occ_pattern
|
||||
! s = 0
|
||||
! do k=1,N_int
|
||||
! if((psi_occ_pattern(k,1,j) /= psi_occ_pattern(k,1,i)).or. &
|
||||
! (psi_occ_pattern(k,2,j) /= psi_occ_pattern(k,2,i))) then
|
||||
! s=1
|
||||
! exit
|
||||
! endif
|
||||
! enddo
|
||||
! if ( s == 0 ) then
|
||||
! print *, 'Error : occ ', j, 'already in wf'
|
||||
! call debug_det(psi_occ_pattern(1,1,j),N_int)
|
||||
! stop
|
||||
! endif
|
||||
! enddo
|
||||
! enddo
|
||||
! !TODO DEBUG
|
||||
|
||||
END_PROVIDER
|
||||
|
||||
subroutine make_s2_eigenfunction
|
||||
@ -253,7 +213,7 @@ subroutine make_s2_eigenfunction
|
||||
integer :: i,j,k
|
||||
integer :: smax, s
|
||||
integer(bit_kind), allocatable :: d(:,:,:), det_buffer(:,:,:)
|
||||
integer :: N_det_new
|
||||
integer :: N_det_new, iproc
|
||||
integer, parameter :: bufsze = 1000
|
||||
logical, external :: is_in_wavefunction
|
||||
|
||||
|
@ -152,7 +152,7 @@ BEGIN_PROVIDER [ double precision, ao_integrals_cache, (0:64*64*64*64) ]
|
||||
END_PROVIDER
|
||||
|
||||
|
||||
double precision function get_ao_bielec_integral(i,j,k,l,map)
|
||||
double precision function get_ao_bielec_integral(i,j,k,l,map) result(result)
|
||||
use map_module
|
||||
implicit none
|
||||
BEGIN_DOC
|
||||
@ -179,15 +179,16 @@ double precision function get_ao_bielec_integral(i,j,k,l,map)
|
||||
call bielec_integrals_index(i,j,k,l,idx)
|
||||
!DIR$ FORCEINLINE
|
||||
call map_get(map,idx,tmp)
|
||||
get_ao_bielec_integral = dble(tmp)
|
||||
tmp = tmp
|
||||
else
|
||||
ii = l-ao_integrals_cache_min
|
||||
ii = ior( ishft(ii,6), k-ao_integrals_cache_min)
|
||||
ii = ior( ishft(ii,6), j-ao_integrals_cache_min)
|
||||
ii = ior( ishft(ii,6), i-ao_integrals_cache_min)
|
||||
get_ao_bielec_integral = ao_integrals_cache(ii)
|
||||
tmp = ao_integrals_cache(ii)
|
||||
endif
|
||||
endif
|
||||
result = tmp
|
||||
end
|
||||
|
||||
|
||||
@ -676,6 +677,7 @@ integer function load_$ao_integrals(filename)
|
||||
real(integral_kind), pointer :: val(:)
|
||||
integer :: iknd, kknd
|
||||
integer*8 :: n, j
|
||||
double precision :: get_$ao_bielec_integral
|
||||
load_$ao_integrals = 1
|
||||
open(unit=66,file=filename,FORM='unformatted',STATUS='UNKNOWN')
|
||||
read(66,err=98,end=98) iknd, kknd
|
||||
|
Loading…
Reference in New Issue
Block a user