mirror of
https://github.com/LCPQ/quantum_package
synced 2024-11-08 07:03:57 +01:00
tri alpha/beta simple pour filter_connected_davidson
This commit is contained in:
parent
b6000270ca
commit
9186be8c56
@ -153,6 +153,7 @@ subroutine davidson_diag_hjj_mrcc(dets_in,u_in,H_jj,energies,dim_in,sze,N_st,Nin
|
||||
! Initialization
|
||||
! ==============
|
||||
|
||||
|
||||
dets_in_sorted(:,:,:) = dets_in(:,:,:)
|
||||
call sort_dets_ab(dets_in_sorted, idx, shortcut, sze, Nint)
|
||||
|
||||
@ -467,57 +468,6 @@ subroutine H_u_0_mrcc(v_0,u_0,H_jj,n,keys_tmp,shortcut,sort_idx,Nint,istate)
|
||||
!$OMP END CRITICAL
|
||||
deallocate(idx,vt)
|
||||
!$OMP END PARALLEL
|
||||
!
|
||||
!
|
||||
!
|
||||
! ASSERT (Nint > 0)
|
||||
! ASSERT (Nint == N_int)
|
||||
! ASSERT (n>0)
|
||||
! PROVIDE ref_bitmask_energy delta_ij
|
||||
! integer, parameter :: block_size = 157
|
||||
! !$OMP PARALLEL DEFAULT(NONE) &
|
||||
! !$OMP PRIVATE(i,hij,j,k,idx,jj,ii,vt) &
|
||||
! !$OMP SHARED(n_det_ref,n_det_non_ref,idx_ref,idx_non_ref,n,H_jj,u_0,keys_tmp,Nint,v_0,istate,delta_ij)
|
||||
! !$OMP DO SCHEDULE(static)
|
||||
! do i=1,n
|
||||
! v_0(i) = H_jj(i) * u_0(i)
|
||||
! enddo
|
||||
! !$OMP END DO
|
||||
! allocate(idx(0:n), vt(n))
|
||||
! Vt = 0.d0
|
||||
! !$OMP DO SCHEDULE(guided)
|
||||
! do i=1,n
|
||||
! idx(0) = i
|
||||
! call filter_connected_davidson(keys_tmp,keys_tmp(1,1,i),Nint,i-1,idx)
|
||||
! do jj=1,idx(0)
|
||||
! j = idx(jj)
|
||||
! ! if ( (dabs(u_0(j)) > 1.d-7).or.((dabs(u_0(i)) > 1.d-7)) ) then
|
||||
! call i_H_j(keys_tmp(1,1,j),keys_tmp(1,1,i),Nint,hij)
|
||||
! hij = hij
|
||||
! vt (i) = vt (i) + hij*u_0(j)
|
||||
! vt (j) = vt (j) + hij*u_0(i)
|
||||
! ! endif
|
||||
! enddo
|
||||
! enddo
|
||||
! !$OMP END DO
|
||||
!
|
||||
! !$OMP DO SCHEDULE(guided)
|
||||
! do ii=1,n_det_ref
|
||||
! i = idx_ref(ii)
|
||||
! do jj = 1, n_det_non_ref
|
||||
! j = idx_non_ref(jj)
|
||||
! vt (i) = vt (i) + delta_ij(ii,jj,istate)*u_0(j)
|
||||
! vt (j) = vt (j) + delta_ij(ii,jj,istate)*u_0(i)
|
||||
! enddo
|
||||
! enddo
|
||||
! !$OMP END DO
|
||||
! !$OMP CRITICAL
|
||||
! do i=1,n
|
||||
! v_0(i) = v_0(i) + vt(i)
|
||||
! enddo
|
||||
! !$OMP END CRITICAL
|
||||
! deallocate(idx,vt)
|
||||
! !$OMP END PARALLEL
|
||||
end
|
||||
|
||||
|
||||
|
@ -324,6 +324,7 @@ subroutine davidson_diag_hjj(dets_in,u_in,H_jj,energies,dim_in,sze,N_st,Nint,iun
|
||||
! ----------------------
|
||||
|
||||
do k=1,N_st
|
||||
! call H_u_0_org(W(1,k,iter),U(1,k,iter),H_jj,sze,dets_in,Nint)
|
||||
call H_u_0(W(1,k,iter),U(1,k,iter),H_jj,sze,dets_in_sorted,shortcut,idx,Nint)
|
||||
enddo
|
||||
|
||||
|
@ -198,16 +198,15 @@ subroutine filter_connected_davidson_warp(key1,warp,key2,Nint,sze,idx)
|
||||
end do
|
||||
endloop = min(warp(2,i_alpha), sze)
|
||||
if(exc_a == 4) then
|
||||
do i_beta=warp(1,i_alpha),endloop
|
||||
beta_loop : do i_beta=warp(1,i_alpha),endloop
|
||||
do ni=1,Nint
|
||||
if(key1(ni,2,i_beta) /= key2(ni,2)) then
|
||||
exit
|
||||
else if(ni == Nint) then
|
||||
idx(l) = i_beta
|
||||
l = l + 1
|
||||
cycle beta_loop
|
||||
end if
|
||||
end do
|
||||
end do
|
||||
idx(l) = i_beta
|
||||
l = l + 1
|
||||
end do beta_loop
|
||||
else
|
||||
do i_beta=warp(1,i_alpha),endloop
|
||||
exc_b = 0
|
||||
@ -283,7 +282,6 @@ subroutine filter_connected_davidson_shortcut(key1,shortcut,key2,Nint,sze,idx)
|
||||
if (Nint==1) then
|
||||
do while(shortcut(i_alpha+1) < sze)
|
||||
i_alpha = i_alpha + 1
|
||||
!print *, i_alpha, shortcut(i_alpha), sze
|
||||
exc_a = popcnt(xor(key1(1,1,shortcut(i_alpha)), key2(1,1)))
|
||||
if(exc_a > 4) then
|
||||
cycle
|
||||
@ -338,6 +336,7 @@ subroutine filter_connected_davidson(key1,key2,Nint,sze,idx)
|
||||
integer*8 :: itmp
|
||||
|
||||
PROVIDE N_con_int det_connections
|
||||
|
||||
ASSERT (Nint > 0)
|
||||
ASSERT (sze >= 0)
|
||||
|
||||
|
@ -138,51 +138,6 @@ subroutine decode_exc(exc,degree,h1,p1,h2,p2,s1,s2)
|
||||
end select
|
||||
end
|
||||
|
||||
!
|
||||
! subroutine get_double_excitation(det1, det2, exc, phase, Nint)
|
||||
! integer, intent(in) :: Nint
|
||||
! integer(bit_kind), intent(in) :: det1(Nint,2)
|
||||
! integer(bit_kind), intent(in) :: det2(Nint,2)
|
||||
! integer, intent(out) :: exc(0:2,2,2)
|
||||
! double precision, intent(out) :: phase
|
||||
!
|
||||
! integer(bit_kind) :: detp, deth, detp_f, deth_f, tmp
|
||||
! integer :: i,sp,hole_fnd
|
||||
!
|
||||
! fi(:) = 0
|
||||
! fs(:) = 0
|
||||
! fn = 0
|
||||
! detp = 0
|
||||
! deth = 0
|
||||
! tmp = 0
|
||||
! deth_f = 0
|
||||
! detp_f = 0
|
||||
! do sp=1,2
|
||||
! do i=1,Nint
|
||||
! tmp = xor(det1(i,sp), det2(i,sp))
|
||||
! if(tmp == 0) then
|
||||
! cycle
|
||||
! end if
|
||||
! detp = iand(det2(i,sp), tmp)
|
||||
! deth = iand(det1(i,sp), tmp)
|
||||
!
|
||||
! detp_f = ior(detp_f, detp)
|
||||
! deth_f = ior(deth_f, deth)
|
||||
!
|
||||
! if(detp /= 0) then
|
||||
! fn = fn + 1
|
||||
! fs(fn) = sp
|
||||
! fi(fn) = i
|
||||
! end if
|
||||
! end do
|
||||
! end do
|
||||
! if(fn /= 2) then
|
||||
! print *, "WHUUUUT??"
|
||||
! stop 1
|
||||
! end if
|
||||
!
|
||||
! end subroutine
|
||||
|
||||
|
||||
subroutine get_double_excitation(det1,det2,exc,phase,Nint)
|
||||
use bitmasks
|
||||
@ -1296,7 +1251,7 @@ subroutine H_u_0(v_0,u_0,H_jj,n,keys_tmp,shortcut,sort_idx,Nint)
|
||||
allocate(idx(0:n), vt(n))
|
||||
Vt = 0.d0
|
||||
v_0 = 0.d0
|
||||
!$OMP DO SCHEDULE(guided)
|
||||
!$OMP DO SCHEDULE(dynamic)
|
||||
|
||||
|
||||
do sh=1,shortcut(0)
|
||||
@ -1344,6 +1299,74 @@ subroutine H_u_0(v_0,u_0,H_jj,n,keys_tmp,shortcut,sort_idx,Nint)
|
||||
end
|
||||
|
||||
|
||||
subroutine H_u_0_org(v_0,u_0,H_jj,n,keys_tmp,Nint)
|
||||
use bitmasks
|
||||
implicit none
|
||||
BEGIN_DOC
|
||||
! Computes v_0 = H|u_0>
|
||||
!
|
||||
! n : number of determinants
|
||||
!
|
||||
! H_jj : array of <j|H|j>
|
||||
END_DOC
|
||||
integer, intent(in) :: n,Nint
|
||||
double precision, intent(out) :: v_0(n)
|
||||
double precision, intent(in) :: u_0(n)
|
||||
double precision, intent(in) :: H_jj(n)
|
||||
integer(bit_kind),intent(in) :: keys_tmp(Nint,2,n)
|
||||
integer, allocatable :: idx(:)
|
||||
double precision :: hij
|
||||
double precision, allocatable :: vt(:)
|
||||
integer :: i,j,k,l, jj,ii,sh
|
||||
integer :: i0, j0
|
||||
|
||||
|
||||
|
||||
ASSERT (Nint > 0)
|
||||
ASSERT (Nint == N_int)
|
||||
ASSERT (n>0)
|
||||
PROVIDE ref_bitmask_energy
|
||||
!$OMP PARALLEL DEFAULT(NONE) &
|
||||
!$OMP PRIVATE(i,hij,j,k,idx,jj,vt,ii) &
|
||||
!$OMP SHARED(n,H_jj,u_0,keys_tmp,Nint,v_0,davidson_threshold)
|
||||
allocate(idx(0:n), vt(n))
|
||||
Vt = 0.d0
|
||||
v_0 = 0.d0
|
||||
!$OMP DO SCHEDULE(guided)
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
do ii=1,n
|
||||
idx(0) = ii
|
||||
i = ii
|
||||
call filter_connected_davidson(keys_tmp,keys_tmp(1,1,ii),Nint,ii-1,idx)
|
||||
|
||||
do jj=1,idx(0)
|
||||
j = idx(jj)
|
||||
if ( dabs(u_0(j)) + dabs(u_0(i)) > davidson_threshold ) then
|
||||
call i_H_j(keys_tmp(1,1,idx(jj)),keys_tmp(1,1,ii),Nint,hij)
|
||||
vt (i) = vt (i) + hij*u_0(j)
|
||||
vt (j) = vt (j) + hij*u_0(i)
|
||||
endif
|
||||
enddo
|
||||
enddo
|
||||
|
||||
!$OMP END DO
|
||||
!$OMP CRITICAL
|
||||
do i=1,n
|
||||
v_0(i) = v_0(i) + vt(i)
|
||||
enddo
|
||||
!$OMP END CRITICAL
|
||||
deallocate(idx,vt)
|
||||
!$OMP END PARALLEL
|
||||
do i=1,n
|
||||
v_0(i) += H_jj(i) * u_0(i)
|
||||
enddo
|
||||
end
|
||||
|
||||
|
||||
|
||||
BEGIN_PROVIDER [ integer, N_con_int ]
|
||||
implicit none
|
||||
|
Loading…
Reference in New Issue
Block a user