1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-07-26 12:47:30 +02:00
qp_plugins_scemama/devel/svdwf/psiSVD_naivBbyB_v1.irp.f

384 lines
12 KiB
FortranFixed
Raw Normal View History

2021-07-28 17:19:18 +02:00
program psiSVD_naivBbyB_v1
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
integer :: i, j, k, l, m, n
double precision :: x, y, h12
double precision, allocatable :: Uref(:,:), Dref(:), Vtref(:,:), Aref(:,:), Vref(:,:)
integer :: rank_max, n_TSVD, n_selected
double precision :: E0, overlop
double precision, allocatable :: H0(:,:)
double precision, allocatable :: eigvec0(:,:), eigval0(:)
integer, allocatable :: numalpha_selected(:), numbeta_selected(:)
integer :: ii
integer :: na_new, nb_new, ind_new, ind_gs
double precision, allocatable :: epsil(:), epsil_energ(:), check_ov(:)
double precision, allocatable :: Uezfio(:,:,:), Dezfio(:,:), Vezfio(:,:,:)
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
integer :: nb_taches
!$OMP PARALLEL
nb_taches = OMP_GET_NUM_THREADS()
!$OMP END PARALLEL
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)
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)
! ---------------------------------------------------------------------------------------
! construct the initial CISD matrix
print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~'
print *, ' CI matrix:', n_det_alpha_unique,'x',n_det_beta_unique
print *, ' N det :', N_det
print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~'
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_alpha_unique) )
allocate( Dref(min(n_det_alpha_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)
allocate( Vref(n_det_beta_unique,n_det_beta_unique) )
do l = 1, n_det_beta_unique
do i = 1, n_det_beta_unique
Vref(i,l) = Vtref(l,i)
enddo
enddo
deallocate( Vtref )
deallocate( Aref )
! ---------------------------------------------------------------------------------------
! ---------------------------------------------------------------------------------------
! initial wavefunction: psi_0
n_TSVD = 1
n_selected = 1
allocate( numalpha_selected(n_selected) , numbeta_selected(n_selected) )
numalpha_selected(1) = 1
numbeta_selected (1) = 1
! get E0 = < psi_0 | H | psi_0 >
allocate( H0(n_selected,n_selected) )
call const_psiHpsi(n_TSVD, n_selected, Uref, Vref, numalpha_selected, numbeta_selected, H0)
E0 = H0(1,1) + nuclear_repulsion
print*, ' ground state E0 = ', E0
deallocate( H0 )
! ---------------------------------------------------------------------------------------
!________________________________________________________________________________________________________
!
! increase the size of psi0 iteratively
!________________________________________________________________________________________________________
rank_max = n_det_alpha_unique*n_det_beta_unique ! 15*15
do while( n_selected .lt. rank_max )
call SYSTEM_CLOCK(COUNT=W_tbeg_it, COUNT_RATE=W_ir)
print*, ' '
print*, ' new iteration '
deallocate( numalpha_selected, numbeta_selected )
n_TSVD = n_TSVD + 1
n_selected = n_TSVD * n_TSVD
allocate( numalpha_selected(n_selected) , numbeta_selected(n_selected) )
l = 0
do i = 1, n_TSVD
do j = 1, n_TSVD
l = l + 1
numalpha_selected(l) = i
numbeta_selected (l) = j
enddo
enddo
if( l.ne.n_selected) then
print *, "error in numbering"
stop
endif
! construct and diagonalise the hamiltonian < psi_selected | H | psi_selected >
allocate( H0(n_selected,n_selected) )
call const_psiHpsi(n_TSVD, n_selected, Uref, Vref, numalpha_selected, numbeta_selected, H0)
allocate( eigvec0(n_selected,n_selected), eigval0(n_selected) )
call lapack_diag(eigval0, eigvec0, H0, n_selected, n_selected)
! get the postsvd ground state
allocate( check_ov(n_selected) )
do l = 1, n_selected
overlop = 0.d0
do i = 1, n_TSVD
ii = i + (i-1)*n_TSVD
overlop = overlop + eigvec0(ii,l) * Dref(i)
enddo
check_ov(l) = dabs(overlop)
enddo
ind_gs = MAXLOC( check_ov, DIM=1 )
overlop = check_ov(ind_gs)
E0 = eigval0(ind_gs)+nuclear_repulsion
print*, ' space dimen = ', n_selected
print*, ' diag energy = ', E0
print*, ' overlop = ', overlop
deallocate( H0 )
deallocate( check_ov, eigval0, eigvec0 )
! ---------------------------------------------------------------------------------------
write(221, *) n_selected, E0
call SYSTEM_CLOCK(COUNT=W_tend_it, COUNT_RATE=W_ir)
W_tot_time_it = real(W_tend_it-W_tbeg_it, kind=8) / real(W_ir, kind=8)
print*, " "
print*, " elapsed time (min) = ", W_tot_time_it/60.d0
end do
!________________________________________________________________________________________________________
!________________________________________________________________________________________________________
deallocate( Dref )
deallocate( Uref, Vref )
! ***************************************************************************************************
! save to ezfion
!allocate( Uezfio(n_det_alpha_unique,rank0,1), Dezfio(rank0,1), Vezfio(n_det_beta_unique,rank0,1) )
!do l = 1, rank0
! Dezfio(l,1) = coeff_psi(l)
! Uezfio(:,l,1) = U0(:,l)
! Vezfio(:,l,1) = V0(:,l)
!enddo
!call ezfio_set_spindeterminants_n_det(N_det)
!call ezfio_set_spindeterminants_n_states(N_states)
!call ezfio_set_spindeterminants_n_det_alpha(n_det_alpha_unique)
!call ezfio_set_spindeterminants_n_det_beta(n_det_beta_unique)
!call ezfio_set_spindeterminants_psi_coef_matrix_rows(psi_bilinear_matrix_rows)
!call ezfio_set_spindeterminants_psi_coef_matrix_columns(psi_bilinear_matrix_columns)
!call ezfio_set_spindeterminants_psi_coef_matrix_values(psi_bilinear_matrix_values)
!call ezfio_set_spindeterminants_n_svd_coefs(rank0)
!call ezfio_set_spindeterminants_psi_svd_alpha(Uezfio)
!call ezfio_set_spindeterminants_psi_svd_beta(Vezfio )
!call ezfio_set_spindeterminants_psi_svd_coefs(Dezfio)
!deallocate( Uezfio, Dezfio, Vezfio )
! ***************************************************************************************************
call SYSTEM_CLOCK(COUNT=W_tend, COUNT_RATE=W_ir)
W_tot_time = real(W_tend - W_tbeg, kind=8) / real(W_ir, kind=8)
print *, ' ___________________________________________________________________'
print *, ' '
print *, " Execution avec ", nb_taches, " threads"
print *, " total elapsed time (min) = ", W_tot_time/60.d0
print *, ' ___________________________________________________________________'
end
subroutine const_psiHpsi(n_TSVD, n_selected, Uref, Vref, numalpha_selected, numbeta_selected, H0)
implicit none
integer, intent(in) :: n_TSVD, n_selected
integer, intent(in) :: numalpha_selected(n_selected), numbeta_selected(n_selected)
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(out) :: H0(n_selected,n_selected)
integer(bit_kind) :: det1(N_int,2)
integer(bit_kind) :: det2(N_int,2)
integer :: degree, na, nb
integer :: i, j, k, l
integer :: iin, jjn, iim, jjm, n, m
double precision :: h12, x
double precision, allocatable :: Utmp(:,:), Vtmp(:,:)
double precision, allocatable :: tmp1(:,:,:,:), tmp2(:,:,:,:)
na = n_det_alpha_unique
nb = n_det_beta_unique
H0(:,:) = 0.d0
allocate( tmp1(nb,nb,n_TSVD,n_TSVD) )
tmp1(:,:,:,:) = 0.d0
allocate( Utmp(n_TSVD,na) )
do i = 1, na
do iin = 1, n_TSVD
Utmp(iin,i) = Uref(i,iin)
enddo
enddo
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE(iin,iim,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,Utmp,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( h12 .eq. 0.d0) cycle
do iin = 1, n_TSVD
x = Utmp(iin,i) * h12
if( x == 0.d0 ) cycle
do iim = 1, n_TSVD
tmp2(j,l,iim,iin) += Utmp(iim,k) * x
enddo
enddo
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP CRITICAL
do iin = 1, n_TSVD
do iim = 1, n_TSVD
do l = 1, nb
do j = 1, nb
tmp1(j,l,iim,iin) += tmp2(j,l,iim,iin)
enddo
enddo
enddo
enddo
!$OMP END CRITICAL
deallocate( tmp2 )
!$OMP END PARALLEL
deallocate( Utmp )
allocate( Vtmp(nb,n_TSVD) )
do iin = 1, n_TSVD
do i = 1, nb
Vtmp(i,iin) = Vref(i,iin)
enddo
enddo
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) &
, Vtmp, size(Vtmp,1) &
, 0.d0, tmp2, size(tmp2,1)*size(tmp2,2)*size(tmp2,3) )
deallocate(tmp1)
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) &
, Vtmp, size(Vtmp,1) &
, 0.d0, tmp1, size(tmp1,1)*size(tmp1,2)*size(tmp1,3) )
deallocate( tmp2, Vtmp )
do n = 1, n_selected
iin = numalpha_selected(n)
jjn = numbeta_selected (n)
do m = 1, n_selected
iim = numalpha_selected(m)
jjm = numbeta_selected (m)
H0(m,n) = tmp1(iin,iim,jjn,jjm)
enddo
enddo
deallocate( tmp1 )
return
end subroutine const_psiHpsi