mirror of
https://gitlab.com/scemama/qp_plugins_scemama.git
synced 2024-11-07 06:33:40 +01:00
335 lines
9.8 KiB
Fortran
335 lines
9.8 KiB
Fortran
program psiSVD_naivBbyB_v2
|
|
|
|
BEGIN_DOC
|
|
! diagonal H in the FSVD space
|
|
END_DOC
|
|
|
|
implicit none
|
|
|
|
read_wf = .True.
|
|
TOUCH read_wf
|
|
|
|
PROVIDE N_int
|
|
|
|
call run()
|
|
end
|
|
|
|
|
|
subroutine run
|
|
|
|
implicit none
|
|
|
|
integer(bit_kind) :: det1(N_int,2), det2(N_int,2)
|
|
integer :: degree, i_state
|
|
double precision :: h12, tol_h
|
|
|
|
integer :: na, nb
|
|
integer :: n_FSVD, n_TSVD
|
|
integer :: i, j, k, l
|
|
|
|
double precision, allocatable :: Uref(:,:), Dref(:), Vtref(:,:), Aref(:,:), Vref(:,:)
|
|
|
|
double precision :: t_beg, t_end
|
|
|
|
call wall_time(t_beg)
|
|
|
|
i_state = 1
|
|
det1(:,1) = psi_det_alpha_unique(:,1)
|
|
det2(:,1) = psi_det_alpha_unique(:,1)
|
|
det1(:,2) = psi_det_beta_unique(:,1)
|
|
det2(:,2) = psi_det_beta_unique(:,1)
|
|
call get_excitation_degree_spin(det1(1,1),det2(1,1),degree,N_int)
|
|
call get_excitation_degree(det1,det2,degree,N_int)
|
|
call i_H_j(det1, det2, N_int, h12)
|
|
|
|
na = n_det_alpha_unique
|
|
nb = n_det_beta_unique
|
|
n_FSVD = min(na,nb)
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
! construct the initial CI matrix
|
|
!
|
|
print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~'
|
|
print *, ' CI matrix:', na, 'x', nb
|
|
print *, ' N det :', N_det
|
|
print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~'
|
|
|
|
allocate( Aref(na,nb) )
|
|
Aref(:,:) = 0.d0
|
|
do k = 1, N_det
|
|
i = psi_bilinear_matrix_rows(k)
|
|
j = psi_bilinear_matrix_columns(k)
|
|
Aref(i,j) = psi_bilinear_matrix_values(k,i_state)
|
|
enddo
|
|
!
|
|
! ---------------------------------------------------------------------------------------
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
! perform a Full SVD
|
|
!
|
|
allocate( Uref(na,na), Dref(n_FSVD), Vtref(nb,nb) )
|
|
call svd_s(Aref, size(Aref,1), Uref, size(Uref,1), Dref, Vtref, size(Vtref,1), na, nb)
|
|
deallocate( Aref )
|
|
|
|
allocate( Vref(n_det_beta_unique,n_det_beta_unique) )
|
|
do l = 1, nb
|
|
do i = 1, nb
|
|
Vref(i,l) = Vtref(l,i)
|
|
enddo
|
|
enddo
|
|
deallocate( Vtref )
|
|
!
|
|
! ---------------------------------------------------------------------------------------
|
|
|
|
|
|
tol_h = 1d-12
|
|
do i = 2, n_FSVD, 10
|
|
call diag_hsvd(i, n_FSVD, tol_h, Uref, Vref, Dref)
|
|
enddo
|
|
call diag_hsvd(n_FSVD, n_FSVD, tol_h, Uref, Vref, Dref)
|
|
|
|
deallocate( Dref, Uref, Vref )
|
|
|
|
call wall_time(t_end)
|
|
print *, ' '
|
|
print *, " total elapsed time (min) = ", (t_end-t_beg)/60.d0
|
|
|
|
end
|
|
|
|
|
|
|
|
! ____________________________________________________________________________________________________________
|
|
!
|
|
! H_svd(a,b,c,d) = sum_{i,j,k,l} H_det(i,j,k,l) U(i,a) V(j,b) U(k,c) V(l,d)
|
|
!
|
|
subroutine diag_hsvd(n_TSVD, n_FSVD, tol_h, Uref, Vref, Dref)
|
|
|
|
implicit none
|
|
|
|
integer, intent(in) :: n_TSVD, n_FSVD
|
|
double precision, intent(in) :: tol_h
|
|
double precision, intent(in) :: Uref(n_det_alpha_unique,n_det_alpha_unique)
|
|
double precision, intent(in) :: Vref(n_det_beta_unique ,n_det_beta_unique)
|
|
double precision, intent(in) :: Dref(n_FSVD)
|
|
|
|
integer(bit_kind) :: det1(N_int,2), det2(N_int,2)
|
|
integer :: degree
|
|
double precision :: h12, x
|
|
|
|
integer :: na, nb, n_selec
|
|
integer :: i, j, k, l
|
|
integer :: a, b, c, d
|
|
|
|
double precision, allocatable :: H_SVD(:,:,:,:)
|
|
double precision, allocatable :: tmp1(:,:,:,:), tmp2(:,:,:,:)
|
|
|
|
integer :: ii, ind_gs
|
|
double precision :: E0, overlop
|
|
double precision, allocatable :: eigvec(:,:), eigval(:), check_ov(:)
|
|
|
|
double precision :: t1, t2, t3
|
|
|
|
na = n_det_alpha_unique
|
|
nb = n_det_beta_unique
|
|
|
|
n_selec = n_TSVD * n_TSVD
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
!
|
|
call wall_time(t1)
|
|
print *, ''
|
|
print *, ' start const psiHpsi in SVD space ...'
|
|
|
|
allocate( tmp1(nb,nb,n_TSVD,n_TSVD) )
|
|
tmp1(:,:,:,:) = 0.d0
|
|
|
|
!$OMP PARALLEL DEFAULT(NONE) &
|
|
!$OMP PRIVATE(a, c, i, j, k, l, h12, x, det1, det2, degree, tmp2) &
|
|
!$OMP SHARED(na, nb, psi_det_alpha_unique, psi_det_beta_unique, &
|
|
!$OMP N_int, n_TSVD, tol_h, Uref, tmp1)
|
|
allocate( tmp2(nb,nb,n_TSVD,n_TSVD) )
|
|
tmp2(:,:,:,:) = 0.d0
|
|
!$OMP DO COLLAPSE(2) SCHEDULE(DYNAMIC,8)
|
|
do l = 1, nb
|
|
do j = 1, nb
|
|
det2(:,2) = psi_det_beta_unique(:,l)
|
|
det1(:,2) = psi_det_beta_unique(:,j)
|
|
call get_excitation_degree_spin(det1(1,2),det2(1,2),degree,N_int)
|
|
if(degree .gt. 2) cycle
|
|
|
|
do k = 1, na
|
|
det2(:,1) = psi_det_alpha_unique(:,k)
|
|
do i = 1, na
|
|
det1(:,1) = psi_det_alpha_unique(:,i)
|
|
call get_excitation_degree(det1,det2,degree,N_int)
|
|
if(degree .gt. 2) cycle
|
|
call i_H_j(det1, det2, N_int, h12)
|
|
if(dabs(h12) .lt. tol_h) cycle
|
|
|
|
! tmp2(j,l,a,c) = sum_{i,k} H(i,j,k,l) U(i,a) U(k,c)
|
|
do c = 1, n_TSVD
|
|
x = Uref(k,c) * h12
|
|
if(dabs(x) .lt. tol_h) cycle
|
|
do a = 1, n_TSVD
|
|
tmp2(j,l,a,c) += Uref(i,a) * x
|
|
enddo
|
|
enddo
|
|
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
|
|
! tmp1 = tmp2
|
|
!$OMP CRITICAL
|
|
do a = 1, n_TSVD
|
|
do c = 1, n_TSVD
|
|
do l = 1, nb
|
|
do j = 1, nb
|
|
tmp1(j,l,a,c) += tmp2(j,l,a,c)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
!$OMP END CRITICAL
|
|
deallocate( tmp2 )
|
|
!$OMP END PARALLEL
|
|
|
|
! tmp2_{l,a,c,b} = sum_{j} tmp1_{j,l,a,c} V_{j,b}
|
|
allocate( tmp2(nb,n_TSVD,n_TSVD,n_TSVD) )
|
|
call DGEMM('T', 'N', nb*n_TSVD*n_TSVD, n_TSVD, nb, 1.d0 &
|
|
, tmp1, size(tmp1,1) &
|
|
, Vref, size(Vref,1) &
|
|
, 0.d0, tmp2, size(tmp2,1)*size(tmp2,2)*size(tmp2,3) )
|
|
deallocate(tmp1)
|
|
|
|
! tmp1_{d,a,c,b} = sum_{l} tmp2_{l,a,c,b} V_{l,d}
|
|
allocate( tmp1(n_TSVD,n_TSVD,n_TSVD,n_TSVD) )
|
|
call DGEMM('T', 'N', n_TSVD*n_TSVD*n_TSVD, n_TSVD, nb, 1.d0 &
|
|
, tmp2, size(tmp2,1) &
|
|
, Vref, size(Vref,1) &
|
|
, 0.d0, tmp1, size(tmp1,1)*size(tmp1,2)*size(tmp1,3) )
|
|
deallocate( tmp2 )
|
|
|
|
allocate( H_SVD(n_TSVD,n_TSVD,n_TSVD,n_TSVD) )
|
|
H_SVD(:,:,:,:) = 0.d0
|
|
do d = 1, n_TSVD
|
|
do c = 1, n_TSVD
|
|
do b = 1, n_TSVD
|
|
do a = 1, n_TSVD
|
|
H_SVD(a,b,c,d) = tmp1(d,a,c,b)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
deallocate( tmp1 )
|
|
|
|
call wall_time(t2)
|
|
print *, ' end const psiHpsi in SVD space after (min)', (t2-t1)/60.d0
|
|
print *, ''
|
|
!
|
|
! ---------------------------------------------------------------------------------------
|
|
|
|
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
!
|
|
print *, ''
|
|
print *, ' diagonalize psiHpsi in SVD space ...'
|
|
|
|
allocate( eigvec(n_selec,n_selec), eigval(n_selec) )
|
|
call lapack_diag(eigval, eigvec, H_SVD, n_selec, n_selec)
|
|
deallocate(H_SVD)
|
|
|
|
! get the postsvd ground state
|
|
allocate( check_ov(n_selec) )
|
|
do l = 1, n_selec
|
|
overlop = 0.d0
|
|
do i = 1, n_TSVD
|
|
ii = i + (i-1) * n_TSVD
|
|
overlop = overlop + eigvec(ii,l) * Dref(i)
|
|
enddo
|
|
check_ov(l) = dabs(overlop)
|
|
enddo
|
|
ind_gs = MAXLOC( check_ov, DIM=1 )
|
|
overlop = check_ov(ind_gs)
|
|
E0 = eigval(ind_gs) + nuclear_repulsion
|
|
|
|
print *, ' - - - - - - - - - - - - '
|
|
print *, ' n_TSVD = ', n_TSVD
|
|
print *, ' space dimen = ', n_selec
|
|
|
|
print *, ' gs num = ', ind_gs
|
|
print *, ' overlop = ', overlop
|
|
|
|
print *, ' diag energy = ', E0
|
|
print *, ' eigen vector: '
|
|
do i = 1, n_selec
|
|
print *, eigvec(i,ind_gs)
|
|
enddo
|
|
print *, ' - - - - - - - - - - - - '
|
|
|
|
write(222, *) n_TSVD, n_selec, E0, overlop
|
|
|
|
deallocate( check_ov, eigval, eigvec )
|
|
|
|
call wall_time(t3)
|
|
print *, ' end diag in SVD space after (min)', (t3-t2)/60.d0
|
|
print *, ''
|
|
!
|
|
! ---------------------------------------------------------------------------------------
|
|
|
|
return
|
|
end subroutine diag_hsvd
|
|
! ____________________________________________________________________________________________________________
|
|
! ____________________________________________________________________________________________________________
|
|
|
|
|
|
|
|
|
|
|
|
|
|
! ! ---------------------------------------------------------------------------------------
|
|
! ! H0 in det basis
|
|
!
|
|
! call wall_time(ti)
|
|
! print *, ''
|
|
! print *, ' start const psiHpsi in det space ...'
|
|
!
|
|
! tol_h = 1d-15
|
|
! allocate( H0_det(na,nb,na,nb) )
|
|
! H0_det(:,:,:,:) = 0.d0
|
|
!
|
|
! !$OMP PARALLEL DEFAULT(NONE) &
|
|
! !$OMP PRIVATE(i, j, k, l, h12, det1, det2, degree) &
|
|
! !$OMP SHARED(na, nb, psi_det_alpha_unique,psi_det_beta_unique, &
|
|
! !$OMP N_int, H0, tol_h)
|
|
! !$OMP DO COLLAPSE(2) SCHEDULE(DYNAMIC,8)
|
|
! do l = 1, nb
|
|
! do j = 1, nb
|
|
! det2(:,2) = psi_det_beta_unique(:,l)
|
|
! det1(:,2) = psi_det_beta_unique(:,j)
|
|
! call get_excitation_degree_spin(det1(1,2),det2(1,2),degree,N_int)
|
|
! if(degree .gt. 2) cycle
|
|
! do k = 1, na
|
|
! det2(:,1) = psi_det_alpha_unique(:,k)
|
|
! do i = 1, na
|
|
! det1(:,1) = psi_det_alpha_unique(:,i)
|
|
! call get_excitation_degree(det1,det2,degree,N_int)
|
|
! if(degree .gt. 2) cycle
|
|
! call i_H_j(det1, det2, N_int, h12)
|
|
! if(dabs(h12) .lt. tol_h) cycle
|
|
! H0(i,j,k,l) = h12
|
|
! enddo
|
|
! enddo
|
|
! enddo
|
|
! enddo
|
|
! !$OMP END DO
|
|
! !$OMP END PARALLEL
|
|
!
|
|
! call wall_time(tf)
|
|
! print *, ' end const psiHpsi in det space after (min)', (ti-tf)/60.d0
|
|
! print *, ''
|
|
! ! ---------------------------------------------------------------------------------------
|