mirror of
https://gitlab.com/scemama/qp_plugins_scemama.git
synced 2024-11-07 22:53:42 +01:00
823 lines
26 KiB
Fortran
823 lines
26 KiB
Fortran
program buildpsi_diagSVDit_v0
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
! perturbative approach to build psi_postsvd
|
|
END_DOC
|
|
|
|
read_wf = .True.
|
|
TOUCH read_wf
|
|
|
|
PROVIDE N_int
|
|
|
|
call run()
|
|
end
|
|
|
|
|
|
subroutine run
|
|
|
|
USE OMP_LIB
|
|
|
|
implicit none
|
|
|
|
integer(bit_kind) :: det1(N_int,2), det2(N_int,2)
|
|
integer :: degree, i_state
|
|
double precision :: h12
|
|
|
|
integer :: i, j, k, l, ii, jj, nn, na, nb
|
|
|
|
double precision :: norm_psi, inv_sqrt_norm_psi
|
|
double precision, allocatable :: Uref(:,:), Dref(:), Vtref(:,:), Aref(:,:), Vref(:,:)
|
|
|
|
double precision :: E0_av, E0_ap, E0pt2
|
|
double precision :: err0, err_tmp, e_tmp, E0, overlop, E0_old, tol_energy
|
|
double precision :: ctmp, htmp, Ept2
|
|
double precision :: E0_postsvd, overlop_postsvd
|
|
double precision :: norm_coeff_psi, inv_sqrt_norm_coeff_psi
|
|
double precision :: overlopU, overlopU_mat, overlopV, overlopV_mat, overlop_psi
|
|
|
|
double precision, allocatable :: Hdiag(:), Hkl(:,:), H0(:,:)
|
|
double precision, allocatable :: psi_postsvd(:), coeff_psi(:), coeff_psi_perturb(:)
|
|
|
|
integer :: n_FSVD, n_selected, n_toselect, n_tmp, it_svd, it_svd_max
|
|
integer :: n_selected2
|
|
integer, allocatable :: numalpha_selected(:), numbeta_selected(:)
|
|
integer, allocatable :: numalpha_toselect(:), numbeta_toselect(:)
|
|
integer, allocatable :: numalpha_tmp(:), numbeta_tmp(:)
|
|
|
|
integer(kind=8) :: W_tbeg, W_tend, W_tbeg_it, W_tend_it, W_ir
|
|
real(kind=8) :: W_tot_time, W_tot_time_it
|
|
real(kind=8) :: CPU_tbeg, CPU_tend, CPU_tbeg_it, CPU_tend_it
|
|
real(kind=8) :: CPU_tot_time, CPU_tot_time_it
|
|
real(kind=8) :: speedup, speedup_it
|
|
integer :: nb_taches
|
|
|
|
!$OMP PARALLEL
|
|
nb_taches = OMP_GET_NUM_THREADS()
|
|
!$OMP END PARALLEL
|
|
|
|
call CPU_TIME(CPU_tbeg)
|
|
call SYSTEM_CLOCK(COUNT=W_tbeg, COUNT_RATE=W_ir)
|
|
|
|
i_state = 1
|
|
|
|
|
|
det1(:,1) = psi_det_alpha_unique(:,1)
|
|
det2(:,1) = psi_det_alpha_unique(:,1)
|
|
call get_excitation_degree_spin(det1(1,1),det2(1,1),degree,N_int)
|
|
det1(:,2) = psi_det_beta_unique(:,1)
|
|
det2(:,2) = psi_det_beta_unique(:,1)
|
|
call get_excitation_degree(det1,det2,degree,N_int)
|
|
call i_H_j(det1, det2, N_int, h12)
|
|
! ---------------------------------------------------------------------------------------
|
|
! construct the initial CISD matrix
|
|
|
|
print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~'
|
|
print *, ' CISD matrix:', n_det_alpha_unique,'x',n_det_beta_unique
|
|
print *, ' N det :', N_det
|
|
print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~'
|
|
|
|
norm_psi = 0.d0
|
|
do k = 1, N_det
|
|
norm_psi = norm_psi + psi_bilinear_matrix_values(k,i_state) &
|
|
* psi_bilinear_matrix_values(k,i_state)
|
|
enddo
|
|
print *, ' initial norm = ', norm_psi
|
|
|
|
allocate( Aref(n_det_alpha_unique,n_det_beta_unique) )
|
|
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(n_det_alpha_unique,n_det_beta_unique) )
|
|
allocate( Dref(n_det_beta_unique) )
|
|
allocate( Vref(n_det_beta_unique,n_det_beta_unique) )
|
|
allocate( Vtref(n_det_beta_unique,n_det_beta_unique) )
|
|
call svd_s(Aref, size(Aref,1), Uref, size(Uref,1), Dref, Vtref, size(Vtref,1) &
|
|
, n_det_alpha_unique, n_det_beta_unique)
|
|
|
|
print *, ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~ '
|
|
print *, ' --- First SVD: ok --- '
|
|
print *, ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~ '
|
|
|
|
do l = 1, n_det_beta_unique
|
|
do i = 1, n_det_beta_unique
|
|
Vref(i,l) = Vtref(l,i)
|
|
enddo
|
|
enddo
|
|
deallocate( Vtref )
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
|
|
nn = n_det_beta_unique
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
! numerote vectors
|
|
|
|
! Full rank
|
|
n_FSVD = nn * nn
|
|
print*, ' Full psi space rank = ', n_FSVD
|
|
|
|
! Truncated rank
|
|
n_selected = 20
|
|
n_selected2 = n_selected * n_selected
|
|
print*, ' initial psi space rank = ', n_selected
|
|
|
|
allocate( numalpha_selected(n_selected) , numbeta_selected(n_selected) )
|
|
do i = 1, n_selected
|
|
numalpha_selected(i) = i
|
|
numbeta_selected (i) = i
|
|
enddo
|
|
! check SVD error
|
|
err0 = 0.d0
|
|
do j = 1, nn
|
|
do i = 1, n_det_alpha_unique
|
|
err_tmp = 0.d0
|
|
do l = 1, n_selected
|
|
err_tmp = err_tmp + Dref(l) * Uref(i,l) * Vref(j,l)
|
|
enddo
|
|
err_tmp = Aref(i,j) - err_tmp
|
|
err0 += err_tmp * err_tmp
|
|
enddo
|
|
enddo
|
|
print *, ' SVD err (%) = ', 100.d0 * dsqrt(err0/norm_psi)
|
|
|
|
deallocate( Aref )
|
|
|
|
! perturbative space rank
|
|
k = 0
|
|
n_toselect = nn*nn - n_selected*n_selected
|
|
allocate( numalpha_toselect(n_toselect) , numbeta_toselect(n_toselect) )
|
|
! nondiagonal blocs
|
|
do i = 1, n_selected
|
|
do j = n_selected+1, nn
|
|
k = k + 1
|
|
numalpha_toselect(k) = j
|
|
numbeta_toselect (k) = i
|
|
enddo
|
|
enddo
|
|
do j = 1, n_selected
|
|
do i = n_selected+1, nn
|
|
k = k + 1
|
|
numalpha_toselect(k) = j
|
|
numbeta_toselect (k) = i
|
|
enddo
|
|
enddo
|
|
! diagonal bloc
|
|
do i = n_selected+1, nn
|
|
do j = n_selected+1, nn
|
|
k = k + 1
|
|
numalpha_toselect(k) = j
|
|
numbeta_toselect (k) = i
|
|
enddo
|
|
enddo
|
|
|
|
if( k.ne.n_toselect ) then
|
|
print*, ' error in numeroting '
|
|
stop
|
|
endif
|
|
print*, ' perturbative psi space rank = ', n_toselect
|
|
|
|
! ---------------------------------------------------------------------------------------
|
|
|
|
|
|
|
|
!________________________________________________________________________________________________________
|
|
!
|
|
! loop over SVD iterations
|
|
!________________________________________________________________________________________________________
|
|
|
|
E0_old = 0.d0
|
|
tol_energy = 1.d0
|
|
it_svd = 0
|
|
it_svd_max = 100
|
|
|
|
do while( ( it_svd .lt. it_svd_max) .and. ( tol_energy .gt. 1d-8 ) )
|
|
|
|
call CPU_TIME(CPU_tbeg_it)
|
|
call SYSTEM_CLOCK(COUNT=W_tbeg_it, COUNT_RATE=W_ir)
|
|
|
|
it_svd = it_svd + 1
|
|
print*, '+ + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + +'
|
|
print*, ' '
|
|
print*, ' iteration', it_svd
|
|
|
|
norm_coeff_psi = 0.d0
|
|
do j = 1, n_selected
|
|
norm_coeff_psi += Dref(j) * Dref(j)
|
|
enddo
|
|
inv_sqrt_norm_coeff_psi = 1.d0 / dsqrt(norm_coeff_psi)
|
|
do j = 1, n_selected
|
|
Dref(j) = Dref(j) * inv_sqrt_norm_coeff_psi
|
|
enddo
|
|
|
|
allocate( H0(n_selected2,n_selected2) )
|
|
print *, ''
|
|
print *, ''
|
|
print *, ''
|
|
print *, '-- Compute H --'
|
|
call const_psihpsi_postsvd_H0(n_selected, n_selected2, Uref, Vref, H0)
|
|
|
|
! avant SVD
|
|
E0 = 0.d0
|
|
do i = 1, n_selected
|
|
ii = (i-1)*n_selected + i
|
|
do j = 1, n_selected
|
|
jj = (j-1)*n_selected + j
|
|
E0 += Dref(j) * H0(jj,ii) * Dref(i)
|
|
enddo
|
|
enddo
|
|
E0_av = E0 + nuclear_repulsion
|
|
print *,' E0 (avant SVD) =', E0_av
|
|
|
|
allocate( psi_postsvd(n_selected2) )
|
|
print *, ' --- Diag post-SVD --- '
|
|
call diag_postsvd(n_selected, n_selected2, Dref, H0, E0_postsvd, overlop_postsvd, psi_postsvd)
|
|
print*, ' postsvd energy = ', E0_postsvd
|
|
deallocate( H0 )
|
|
|
|
! post-SVD
|
|
print *, ' --- SVD --- '
|
|
!Dref(:) = 0.d0
|
|
call perform_newpostSVD(n_selected, n_selected2, psi_postsvd, Uref, Vref, Dref)
|
|
deallocate( psi_postsvd )
|
|
|
|
allocate( H0(n_selected2,n_selected2) )
|
|
print *, ' --- Compute H --- '
|
|
call const_psihpsi_postsvd_H0(n_selected, n_selected2, Uref, Vref, H0)
|
|
|
|
E0 = 0.d0
|
|
norm_coeff_psi = 0.d0
|
|
do i = 1, n_selected
|
|
ii = (i-1)*n_selected + i
|
|
do j = 1, n_selected
|
|
jj = (j-1)*n_selected + j
|
|
E0 += Dref(j) * H0(jj,ii) * Dref(i)
|
|
enddo
|
|
norm_coeff_psi += Dref(i) * Dref(i)
|
|
enddo
|
|
E0_ap = E0 + nuclear_repulsion
|
|
print *,' E0 (apres SVD) =', E0_ap
|
|
!print *,' norm =', norm_coeff_psi
|
|
|
|
deallocate(H0)
|
|
|
|
print *, ' --- Perturbation --- '
|
|
|
|
allocate( Hdiag(n_toselect), Hkl(n_selected2,n_toselect) )
|
|
call const_Hdiag_Hkl(n_selected, n_selected2, n_toselect, Uref, Vref, numalpha_toselect, numbeta_toselect, Hdiag, Hkl)
|
|
!do l = 1, n_toselect
|
|
! na = numalpha_toselect(l)
|
|
! nb = numbeta_toselect (l)
|
|
! print *, na, nb, Hdiag(l)
|
|
!enddo
|
|
|
|
! evaluate the coefficients for all the vectors
|
|
allocate( coeff_psi_perturb(n_toselect) )
|
|
ept2 = 0.d0
|
|
do ii = 1, n_toselect
|
|
!na = numalpha_toselect(ii)
|
|
!nb = numbeta_toselect (ii)
|
|
ctmp = 0.d0
|
|
do i = 1, n_selected
|
|
l = (i-1)*n_selected + i
|
|
ctmp += Dref(i) * Hkl(l,ii)
|
|
enddo
|
|
coeff_psi_perturb(ii) = ctmp / ( E0_ap - (Hdiag(ii)+nuclear_repulsion) )
|
|
ept2 += ctmp*ctmp / ( E0_ap - (Hdiag(ii)+nuclear_repulsion) )
|
|
enddo
|
|
E0pt2 = E0_ap + ept2
|
|
print *, ' perturb energy = ', E0pt2, ept2
|
|
|
|
tol_energy = 100.d0 * dabs(E0pt2-E0_old) / dabs(E0pt2)
|
|
E0_old = E0pt2
|
|
|
|
deallocate( Hdiag, Hkl)
|
|
|
|
print *, ' --- SVD --- '
|
|
call perform_newSVD(n_selected, n_selected2, n_toselect, numalpha_toselect, numbeta_toselect, coeff_psi_perturb, Uref, Vref, Dref)
|
|
|
|
deallocate( coeff_psi_perturb )
|
|
|
|
write(11,'(i5,4x,4(f22.15,2x))') it_svd, E0_av, E0_postsvd, E0_ap, E0pt2
|
|
|
|
call CPU_TIME(CPU_tend_it)
|
|
call SYSTEM_CLOCK(COUNT=W_tend_it, COUNT_RATE=W_ir)
|
|
CPU_tot_time_it = CPU_tend_it - CPU_tbeg_it
|
|
W_tot_time_it = real(W_tend_it-W_tbeg_it, kind=8) / real(W_ir, kind=8)
|
|
speedup_it = CPU_tot_time_it / W_tot_time_it
|
|
print '(//, 3X, "elapsed time = ", 1PE10.3, " min.", /, &
|
|
& 3X, "CPU time = ", 1PE10.3, " min.", /, &
|
|
& 3X, "speed up = ", 1PE10.3,//)', W_tot_time_it/60.d0, CPU_tot_time_it/60.d0, speedup_it
|
|
|
|
end do
|
|
!________________________________________________________________________________________________________
|
|
!________________________________________________________________________________________________________
|
|
|
|
|
|
|
|
deallocate( Uref, Vref, Dref )
|
|
|
|
|
|
call CPU_TIME(CPU_tend)
|
|
call SYSTEM_CLOCK(COUNT=W_tend, COUNT_RATE=W_ir)
|
|
CPU_tot_time = CPU_tend - CPU_tbeg
|
|
W_tot_time = real(W_tend - W_tbeg, kind=8) / real(W_ir, kind=8)
|
|
speedup = CPU_tot_time / W_tot_time
|
|
print *,' ___________________________________________________________________'
|
|
print '(//,3X,"Execution avec ",i2," threads")',nb_taches
|
|
print '(//, 3X, "elapsed time = ", 1PE10.3, " min.", /, &
|
|
& 3X, "CPU time = ", 1PE10.3, " min.", /, &
|
|
& 3X, "speed up = ", 1PE10.3 ,// )', W_tot_time/60.d0, CPU_tot_time/60.d0, speedup
|
|
print *,' ___________________________________________________________________'
|
|
|
|
|
|
end
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
subroutine const_psihpsi_postsvd_H0(n_selected, n_selected2, Uref, Vref, H0)
|
|
|
|
USE OMP_LIB
|
|
|
|
implicit none
|
|
|
|
integer, intent(in) :: n_selected, n_selected2
|
|
double precision, intent(in) :: Uref(n_det_alpha_unique,n_det_beta_unique)
|
|
double precision, intent(in) :: Vref(n_det_beta_unique ,n_det_beta_unique)
|
|
double precision, intent(out) :: H0(n_selected2,n_selected2)
|
|
|
|
integer(bit_kind) :: det1(N_int,2), det2(N_int,2)
|
|
integer :: i, j, k, l, degree
|
|
integer :: ii0, jj0, ii, jj, n, m, np, mp
|
|
integer :: nn0, mm0, nn, mm, ind_gs
|
|
double precision :: h12, x
|
|
|
|
double precision, allocatable :: H0_tmp(:,:)
|
|
|
|
|
|
H0(:,:) = 0.d0
|
|
|
|
!$OMP PARALLEL DEFAULT(NONE) &
|
|
!$OMP PRIVATE(n,np,nn0,nn,ii0,jj0,x,m,mp,mm0,mm,ii,jj,i,j,k,l,h12,det1,det2,H0_tmp,degree) &
|
|
!$OMP SHARED(n_det_alpha_unique,n_det_beta_unique,psi_det_alpha_unique,psi_det_beta_unique, &
|
|
!$OMP N_int,n_selected,n_selected2,Uref,Vref,H0 )
|
|
allocate( H0_tmp(n_selected2,n_selected2) )
|
|
H0_tmp(:,:) = 0.d0
|
|
!$OMP DO COLLAPSE(2) SCHEDULE(DYNAMIC,8)
|
|
do i = 1, n_det_alpha_unique
|
|
do k = 1, n_det_alpha_unique
|
|
det1(:,1) = psi_det_alpha_unique(:,i)
|
|
det2(:,1) = psi_det_alpha_unique(:,k)
|
|
call get_excitation_degree_spin(det1(1,1),det2(1,1),degree,N_int)
|
|
if (degree .gt. 2) then
|
|
cycle
|
|
endif
|
|
do j = 1, n_det_beta_unique
|
|
det1(:,2) = psi_det_beta_unique(:,j)
|
|
do l = 1, n_det_beta_unique
|
|
det2(:,2) = psi_det_beta_unique(:,l)
|
|
call get_excitation_degree(det1,det2,degree,N_int)
|
|
if (degree .gt. 2) then
|
|
cycle
|
|
endif
|
|
! !!!
|
|
call i_H_j(det1, det2, N_int, h12)
|
|
! !!!
|
|
! ~~~ H0 ~~~
|
|
do n = 1, n_selected
|
|
nn0 = (n-1)*n_selected
|
|
do np = 1, n_selected
|
|
nn = nn0 + np
|
|
x = Uref(k,n) * Vref(l,np) * h12
|
|
do m = 1, n_selected
|
|
mm0 = (m-1)*n_selected
|
|
do mp = 1, n_selected
|
|
mm = mm0 + mp
|
|
H0_tmp(mm,nn) += Uref(i,m) * Vref(j,mp) * x
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
! ~~~ ~~~~~~ ~~~
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP CRITICAL
|
|
do n = 1, n_selected2
|
|
do m = 1, n_selected2
|
|
H0(m,n) += H0_tmp(m,n)
|
|
enddo
|
|
enddo
|
|
!$OMP END CRITICAL
|
|
deallocate( H0_tmp )
|
|
!$OMP END PARALLEL
|
|
|
|
return
|
|
end subroutine const_psihpsi_postsvd_H0
|
|
|
|
|
|
|
|
|
|
|
|
subroutine diag_postsvd(n_selected, n_selected2, Dref, H0, E0, overlop, psi_postsvd )
|
|
|
|
USE OMP_LIB
|
|
|
|
implicit none
|
|
|
|
integer, intent(in) :: n_selected, n_selected2
|
|
double precision, intent(in) :: H0(n_selected2,n_selected2)
|
|
double precision, intent(in) :: Dref(n_det_beta_unique)
|
|
double precision, intent(out) :: E0, overlop, psi_postsvd(n_selected2)
|
|
|
|
integer(bit_kind) :: det1(N_int,2), det2(N_int,2)
|
|
integer :: i, j, k, l, degree
|
|
integer :: ii0, jj0, ii, jj, n, m, np, mp
|
|
integer :: nn0, mm0, nn, mm, ind_gs
|
|
double precision :: h12, x
|
|
|
|
double precision, allocatable :: eigvec0(:,:), eigval0(:), check_ov(:)
|
|
|
|
! diagonalize H0
|
|
allocate( eigvec0(n_selected2,n_selected2), eigval0(n_selected2) )
|
|
call lapack_diag(eigval0, eigvec0, H0, n_selected2, n_selected2)
|
|
|
|
! get the postsvd ground state
|
|
allocate( check_ov(n_selected2) )
|
|
do l = 1, n_selected2
|
|
overlop = 0.d0
|
|
do i = 1, n_selected
|
|
ii = n_selected*(i-1) + i
|
|
overlop = overlop + eigvec0(ii,l) * Dref(i)
|
|
enddo
|
|
check_ov(l) = dabs(overlop)
|
|
enddo
|
|
ind_gs = MAXLOC( check_ov, DIM=1 )
|
|
!ind_gs = 1
|
|
overlop = check_ov(ind_gs)
|
|
E0 = eigval0(ind_gs)+nuclear_repulsion
|
|
psi_postsvd = eigvec0(:,ind_gs)
|
|
|
|
deallocate( check_ov, eigvec0, eigval0 )
|
|
|
|
return
|
|
end subroutine diag_postsvd
|
|
|
|
|
|
|
|
|
|
subroutine const_Hdiag_Hkl(n_selected, n_selected2, n_toselect, Uref, Vref, numalpha_toselect, numbeta_toselect, Hdiag, Hkl)
|
|
|
|
implicit none
|
|
|
|
integer, intent(in) :: n_selected, n_selected2, n_toselect
|
|
integer, intent(in) :: numalpha_toselect(n_toselect), numbeta_toselect(n_toselect)
|
|
double precision, intent(in) :: Uref(n_det_alpha_unique,n_det_beta_unique)
|
|
double precision, intent(in) :: Vref(n_det_beta_unique ,n_det_beta_unique)
|
|
double precision, intent(out) :: Hdiag(n_toselect), Hkl(n_selected2,n_toselect)
|
|
|
|
integer(bit_kind) :: det1(N_int,2)
|
|
integer(bit_kind) :: det2(N_int,2)
|
|
integer :: degree
|
|
|
|
integer :: i, j, k, l
|
|
integer :: ii0, jj0, ii, jj, n, m, np, mp
|
|
double precision :: h12, y
|
|
|
|
double precision, allocatable :: Hdiag_tmp(:), Hkl_tmp(:,:)
|
|
|
|
|
|
Hdiag(:) = 0.d0
|
|
Hkl(:,:) = 0.d0
|
|
|
|
!$OMP PARALLEL DEFAULT(NONE) &
|
|
!$OMP PRIVATE(n,ii0,jj0,y,m,mp,ii,jj,i,j,k,l,h12,det1,det2,Hdiag_tmp,Hkl_tmp,degree) &
|
|
!$OMP SHARED(n_det_alpha_unique,n_det_beta_unique,psi_det_alpha_unique,psi_det_beta_unique, &
|
|
!$OMP N_int,n_selected,n_toselect,Uref,Vref,numalpha_toselect,numbeta_toselect, &
|
|
!$OMP Hkl,Hdiag,n_selected2 )
|
|
allocate( Hdiag_tmp(n_toselect), Hkl_tmp(n_selected2,n_toselect) )
|
|
Hdiag_tmp(:) = 0.d0
|
|
Hkl_tmp(:,:) = 0.d0
|
|
!$OMP DO COLLAPSE(2) SCHEDULE(DYNAMIC,8)
|
|
do i = 1, n_det_alpha_unique
|
|
do k = 1, n_det_alpha_unique
|
|
det1(:,1) = psi_det_alpha_unique(:,i)
|
|
det2(:,1) = psi_det_alpha_unique(:,k)
|
|
! !!!
|
|
! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
|
|
call get_excitation_degree_spin(det1(1,1),det2(1,1),degree,N_int)
|
|
if (degree .gt. 2) then
|
|
cycle
|
|
endif
|
|
! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
|
|
! !!!
|
|
do j = 1, n_det_beta_unique
|
|
det1(:,2) = psi_det_beta_unique(:,j)
|
|
do l = 1, n_det_beta_unique
|
|
det2(:,2) = psi_det_beta_unique(:,l)
|
|
! !!!
|
|
! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
|
|
call get_excitation_degree(det1,det2,degree,N_int)
|
|
if (degree .gt. 2) then
|
|
cycle
|
|
endif
|
|
! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
|
|
! !!!
|
|
call i_H_j(det1, det2, N_int, h12)
|
|
! ~ ~ ~ H ~ ~ ~
|
|
do n = 1, n_toselect
|
|
ii0 = numalpha_toselect(n)
|
|
jj0 = numbeta_toselect (n)
|
|
y = Uref(k,ii0) * Vref(l,jj0) * h12
|
|
! Hdiag
|
|
Hdiag_tmp(n) += y * Uref(i,ii0) * Vref(j,jj0)
|
|
do m = 1, n_selected
|
|
ii = (m-1)*n_selected
|
|
do mp = 1, n_selected
|
|
jj = ii + mp
|
|
! Hkl
|
|
Hkl_tmp(jj,n) += Uref(i,m) * Vref(j,mp) * y
|
|
enddo
|
|
enddo
|
|
enddo
|
|
! ~ ~ ~ ! ! ! ~ ~ ~
|
|
enddo
|
|
enddo
|
|
! !!!
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP CRITICAL
|
|
do n = 1, n_toselect
|
|
Hdiag(n) += Hdiag_tmp(n)
|
|
do m = 1, n_selected2
|
|
Hkl(m,n) += Hkl_tmp(m,n)
|
|
enddo
|
|
enddo
|
|
!$OMP END CRITICAL
|
|
deallocate( Hdiag_tmp,Hkl_tmp )
|
|
!$OMP END PARALLEL
|
|
|
|
end subroutine const_Hdiag_Hkl
|
|
|
|
|
|
|
|
|
|
|
|
|
|
subroutine perform_newSVD(n_selected, n_selected2, n_toselect, numalpha_toselect, numbeta_toselect, coeff_psi_perturb, Uref, Vref, Dref)
|
|
|
|
USE OMP_LIB
|
|
|
|
integer, intent(in) :: n_selected, n_toselect, n_selected2
|
|
integer, intent(in) :: numalpha_toselect(n_toselect), numbeta_toselect(n_toselect)
|
|
double precision, intent(in) :: coeff_psi_perturb(n_toselect)
|
|
double precision, intent(inout) :: Uref(n_det_alpha_unique,n_det_beta_unique)
|
|
double precision, intent(inout) :: Vref(n_det_beta_unique ,n_det_beta_unique)
|
|
double precision, intent(inout) :: Dref(n_det_beta_unique)
|
|
|
|
integer :: mm, nn, i, j, ii0, ii, l, jj, na, nb
|
|
double precision :: err0, err_norm, err_tmp, norm_tmp
|
|
double precision :: overlopU_mat, overlopV_mat, overlopU, overlopV
|
|
double precision, allocatable :: S_mat(:,:), SxVt(:,:)
|
|
double precision, allocatable :: U_svd(:,:), V_svd(:,:)
|
|
double precision, allocatable :: U_newsvd(:,:), V_newsvd(:,:), Vt_newsvd(:,:), D_newsvd(:), A_newsvd(:,:)
|
|
|
|
mm = n_det_alpha_unique
|
|
nn = n_det_beta_unique
|
|
|
|
allocate( U_svd(mm,nn) , V_svd(nn,nn) , S_mat(nn,nn) )
|
|
|
|
U_svd(:,:) = Uref(:,:)
|
|
V_svd(:,:) = Vref(:,:)
|
|
S_mat(:,:) = 0.d0
|
|
norm_tmp = 0.d0
|
|
do j = 1, n_det_beta_unique
|
|
S_mat(j,j) = Dref(j)
|
|
norm_tmp += S_mat(j,j)*S_mat(j,j)
|
|
enddo
|
|
do l = 1, n_toselect
|
|
na = numalpha_toselect(l)
|
|
nb = numbeta_toselect (l)
|
|
S_mat(na,nb) = coeff_psi_perturb(l)
|
|
norm_tmp += S_mat(na,nb)*S_mat(na,nb)
|
|
enddo
|
|
|
|
print*, ' norm de S_mat =', norm_tmp
|
|
!norm_tmp = 1.d0/dsqrt(norm_tmp)
|
|
!do i = 1, nn
|
|
! do j = 1, nn
|
|
! S_mat(j,i) = S_mat(j,i) * norm_tmp
|
|
! enddo
|
|
!enddo
|
|
|
|
|
|
! first compute S_mat x transpose(V_svd)
|
|
allocate( SxVt(nn,nn) )
|
|
call dgemm( 'N', 'T', nn, nn, nn, 1.d0 &
|
|
, S_mat , size(S_mat,1) &
|
|
, V_svd , size(V_svd,1) &
|
|
, 0.d0, SxVt, size(SxVt ,1) )
|
|
! then compute U_svd x SxVt
|
|
allocate( A_newsvd(mm,nn) )
|
|
call dgemm( 'N', 'N', mm, nn, nn, 1.d0 &
|
|
, U_svd , size(U_svd ,1) &
|
|
, SxVt , size(SxVt ,1) &
|
|
, 0.d0, A_newsvd, size(A_newsvd,1) )
|
|
deallocate( SxVt )
|
|
|
|
! perform new SVD
|
|
allocate( U_newsvd(mm,nn), Vt_newsvd(nn,nn), D_newsvd(nn) )
|
|
call svd_s( A_newsvd, size(A_newsvd,1), U_newsvd, size(U_newsvd,1), D_newsvd, Vt_newsvd, size(Vt_newsvd,1), mm, nn)
|
|
print *, ' +++ new perturbative SVD is performed +++ '
|
|
allocate( V_newsvd(nn,nn) )
|
|
do l = 1, nn
|
|
do j = 1, nn
|
|
V_newsvd(j,l) = Vt_newsvd(l,j)
|
|
enddo
|
|
enddo
|
|
|
|
! check SVD error
|
|
err0 = 0.d0
|
|
err_norm = 0.d0
|
|
do j = 1, nn
|
|
do i = 1, mm
|
|
err_tmp = 0.d0
|
|
do l = 1, nn
|
|
err_tmp = err_tmp + D_newsvd(l) * U_newsvd(i,l) * V_newsvd(j,l)
|
|
enddo
|
|
err_tmp = A_newsvd(i,j) - err_tmp
|
|
err0 += err_tmp * err_tmp
|
|
err_norm += A_newsvd(i,j) * A_newsvd(i,j)
|
|
enddo
|
|
enddo
|
|
print *, ' SVD err (%) = ', 100.d0 * dsqrt(err0/err_norm)
|
|
print *, ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ '
|
|
|
|
|
|
do l = 1, nn
|
|
Dref(l) = D_newsvd(l)
|
|
Uref(:,l) = U_newsvd(:,l)
|
|
Vref(:,l) = V_newsvd(:,l)
|
|
enddo
|
|
!print *, Dref(:)
|
|
|
|
|
|
overlopU_mat = 0.d0
|
|
overlopV_mat = 0.d0
|
|
do i = 1, nn
|
|
do j = 1, nn
|
|
overlopU = 0.d0
|
|
do ii = 1, mm
|
|
overlopU += Uref(ii,j) * Uref(ii,i)
|
|
enddo
|
|
overlopU_mat += overlopU
|
|
overlopV = 0.d0
|
|
do ii = 1, nn
|
|
overlopV += Vref(ii,j) * Vref(ii,i)
|
|
enddo
|
|
overlopV_mat += overlopV
|
|
enddo
|
|
enddo
|
|
print *, 'overlop U =', overlopU_mat
|
|
print *, 'overlop V =', overlopV_mat
|
|
|
|
|
|
deallocate( U_newsvd, V_newsvd, Vt_newsvd, D_newsvd, A_newsvd )
|
|
|
|
return
|
|
|
|
end subroutine perform_newSVD
|
|
|
|
|
|
|
|
|
|
|
|
subroutine perform_newpostSVD(n_selected, n_selected2, psi_postsvd, Uref, Vref, Dref)
|
|
|
|
USE OMP_LIB
|
|
|
|
integer, intent(in) :: n_selected, n_selected2
|
|
double precision, intent(in) :: psi_postsvd(n_selected2)
|
|
double precision, intent(inout) :: Uref(n_det_alpha_unique,n_det_beta_unique)
|
|
double precision, intent(inout) :: Vref(n_det_beta_unique ,n_det_beta_unique)
|
|
double precision, intent(inout) :: Dref(n_det_beta_unique)
|
|
|
|
integer :: mm, nn, i, j, ii0, ii, l, jj, na, nb
|
|
double precision :: err0, err_norm, err_tmp, norm_tmp
|
|
double precision :: overlopU_mat, overlopV_mat, overlopU, overlopV
|
|
double precision, allocatable :: S_mat(:,:), SxVt(:,:)
|
|
double precision, allocatable :: U_svd(:,:), V_svd(:,:)
|
|
double precision, allocatable :: U_newsvd(:,:), V_newsvd(:,:), Vt_newsvd(:,:), D_newsvd(:), A_newsvd(:,:)
|
|
|
|
mm = n_det_alpha_unique
|
|
nn = n_det_beta_unique
|
|
|
|
allocate( U_svd(mm,n_selected) , V_svd(nn,n_selected) , S_mat(n_selected,n_selected) )
|
|
|
|
U_svd(:,:) = Uref(:,1:n_selected)
|
|
V_svd(:,:) = Vref(:,1:n_selected)
|
|
S_mat(:,:) = 0.d0
|
|
do i = 1, n_selected
|
|
ii = (i-1)*n_selected
|
|
do j = 1, n_selected
|
|
jj = ii + j
|
|
S_mat(i,j) = psi_postsvd(jj)
|
|
enddo
|
|
enddo
|
|
|
|
! first compute S_mat x transpose(V_svd)
|
|
allocate( SxVt(n_selected,nn) )
|
|
call dgemm( 'N', 'T', n_selected, nn, n_selected, 1.d0 &
|
|
, S_mat , size(S_mat,1) &
|
|
, V_svd , size(V_svd,1) &
|
|
, 0.d0, SxVt, size(SxVt ,1) )
|
|
! then compute U_svd x SxVt
|
|
allocate( A_newsvd(mm,nn) )
|
|
call dgemm( 'N', 'N', mm, nn, n_selected, 1.d0 &
|
|
, U_svd , size(U_svd ,1) &
|
|
, SxVt , size(SxVt ,1) &
|
|
, 0.d0, A_newsvd, size(A_newsvd,1) )
|
|
deallocate( SxVt )
|
|
|
|
! perform new SVD
|
|
allocate( U_newsvd(mm,nn), Vt_newsvd(nn,nn), D_newsvd(nn) )
|
|
call svd_s( A_newsvd, size(A_newsvd,1), U_newsvd, size(U_newsvd,1), D_newsvd, Vt_newsvd, size(Vt_newsvd,1), mm, nn)
|
|
print *, ' +++ new SVD is performed +++ '
|
|
allocate( V_newsvd(nn,nn) )
|
|
do l = 1, nn
|
|
do j = 1, nn
|
|
V_newsvd(j,l) = Vt_newsvd(l,j)
|
|
enddo
|
|
enddo
|
|
|
|
! check SVD error
|
|
err0 = 0.d0
|
|
err_norm = 0.d0
|
|
do j = 1, nn
|
|
do i = 1, mm
|
|
err_tmp = 0.d0
|
|
do l = 1, n_selected
|
|
err_tmp = err_tmp + D_newsvd(l) * U_newsvd(i,l) * V_newsvd(j,l)
|
|
enddo
|
|
err_tmp = A_newsvd(i,j) - err_tmp
|
|
err0 += err_tmp * err_tmp
|
|
err_norm += A_newsvd(i,j) * A_newsvd(i,j)
|
|
enddo
|
|
enddo
|
|
print *, ' SVD err (%) = ', 100.d0 * dsqrt(err0/err_norm)
|
|
print *, ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ '
|
|
|
|
do l = 1, n_selected
|
|
Dref(l) = D_newsvd(l)
|
|
Uref(:,l) = U_newsvd(:,l)
|
|
Vref(:,l) = V_newsvd(:,l)
|
|
enddo
|
|
! print *, Dref(:)
|
|
|
|
overlopU_mat = 0.d0
|
|
overlopV_mat = 0.d0
|
|
do i = 1, nn
|
|
do j = 1, nn
|
|
overlopU = 0.d0
|
|
do ii = 1, mm
|
|
overlopU += Uref(ii,j) * Uref(ii,i)
|
|
enddo
|
|
overlopU_mat += overlopU
|
|
overlopV = 0.d0
|
|
do ii = 1, nn
|
|
overlopV += Vref(ii,j) * Vref(ii,i)
|
|
enddo
|
|
overlopV_mat += overlopV
|
|
enddo
|
|
enddo
|
|
print *, 'overlop U =', overlopU_mat
|
|
print *, 'overlop V =', overlopV_mat
|
|
|
|
|
|
deallocate( U_newsvd, V_newsvd, Vt_newsvd, D_newsvd, A_newsvd )
|
|
|
|
return
|
|
|
|
end subroutine perform_newpostSVD
|
|
|
|
|