1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-07-26 04:37:31 +02:00
qp_plugins_scemama/devel/svdwf/SQ_klHkl_v1.irp.f
2021-11-02 16:18:07 +01:00

512 lines
18 KiB
Fortran

program SQ_klHkl_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
USE bitmasks
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
double precision :: E0, overlop, Ept2
double precision, allocatable :: H0(:,:)
double precision, allocatable :: eigvec0(:,:), eigval0(:), coeff_psi(:), coeff_tmp(:)
integer :: ii, ia, ib
double precision, allocatable :: Hdiag(:), Hkl_save(:,:), Hkl_1d(:), Hkl_tmp(:,:), Hdiag_tmp(:)
integer :: na_new, nb_new, ind_new, ind_gs
double precision :: ctmp, coeff_new
double precision, allocatable :: epsil(:), epsil_energ(:), check_ov(:)
double precision, allocatable :: Uezfio(:,:,:), Dezfio(:,:), Vezfio(:,:,:)
integer :: ibeg_alpha, ibeg_beta, iend_alpha, iend_beta
integer :: n_toselect, na_max, nb_max
integer, allocatable :: numalpha_toselect(:), numbeta_toselect(:)
double precision :: t_beg, t_end, ti, tf
integer :: nb_taches
!$OMP PARALLEL
nb_taches = OMP_GET_NUM_THREADS()
!$OMP END PARALLEL
call wall_time(ti)
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 *, ' CISD 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 wall_time(t_beg)
call svd_s(Aref, size(Aref,1), Uref, size(Uref,1), Dref, Vtref &
, size(Vtref,1), n_det_alpha_unique, n_det_beta_unique)
call wall_time(t_end)
print *, " SVD is performed after (min)", (t_end-t_beg)/60.
deallocate( Aref , Dref )
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 )
!
! ---------------------------------------------------------------------------------------
! ---------------------------------------------------------------------------------------
! numerote | k l > toselect
ibeg_alpha = 1
iend_alpha = n_det_alpha_unique
na_max = iend_alpha - ibeg_alpha + 1
ibeg_beta = 1
iend_beta = n_det_beta_unique
nb_max = iend_beta - ibeg_beta + 1
n_toselect = na_max * nb_max
print *, ' na_max = ', na_max
print *, ' nb_max = ', nb_max
print *, ' n_toselect = ', n_toselect
allocate( numalpha_toselect(n_toselect) , numbeta_toselect(n_toselect) )
k = 0
do i = ibeg_alpha, iend_alpha
do j = ibeg_beta, iend_beta
k = k + 1
numalpha_toselect(k) = i
numbeta_toselect (k) = j
enddo
enddo
if( k.ne.n_toselect ) then
print *, " error in numbering"
stop
endif
!
! ---------------------------------------------------------------------------------------
double precision, allocatable :: A_ik(:,:), B_jl(:,:), T(:,:,:,:)
allocate( A_ik(n_det_alpha_unique,n_det_alpha_unique) , B_jl(n_det_beta_unique,n_det_beta_unique) )
allocate( T(n_det_alpha_unique,n_det_alpha_unique,n_det_beta_unique,n_det_beta_unique ) )
call const_AB(A_ik, B_jl)
call const_2b(T)
allocate( Hdiag(n_toselect) )
call const_Hdiag(n_toselect, Uref, Vref, numalpha_toselect, numbeta_toselect, A_ik, B_jl, T, Hdiag)
deallocate( A_ik, B_jl, T )
open(UNIT=11, FILE="SQ_klHkl_v1.dat", ACTION="WRITE")
do i = 1, n_toselect
write(11, '(2(I5,2X), 5X, E15.7)') numalpha_toselect(i), numbeta_toselect(i), Hdiag(i)
enddo
close(11)
deallocate( Hdiag )
deallocate( Uref, Vref )
deallocate( numalpha_toselect, numbeta_toselect )
call wall_time(tf)
print *, ' ___________________________________________________________________'
print *, ' '
!print *, " Execution avec ", nb_taches, " threads"
print *, " Execution avec 1 threads"
print *, " total elapsed time (min) = ", (tf-ti)/60.d0
print *, ' ___________________________________________________________________'
end
!/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\!
!---! !---! !---! !---! !---! !---! !---! !---! !---!
! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\!
!---! !---! !---! !---! !---! !---! !---! !---! !---!
! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\!
!---! !---! !---! !---! !---! !---! !---! !---! !---!
! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\!
!---! !---! !---! !---! !---! !---! !---! !---! !---!
! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !
!___! !___! !___! !___! !___! !___! !___! !___! !___!
!___! !___! !___! !___! !___! !___! !___! !___! !___!
subroutine const_AB(A_ik, B_jl)
USE bitmasks
implicit none
! see one_e_dm_mo_alpha & one_e_dm_mo_beta in determinants/density_matrix.irp.f
double precision, intent(out) :: A_ik(n_det_alpha_unique,n_det_alpha_unique)
double precision, intent(out) :: B_jl(n_det_beta_unique ,n_det_beta_unique )
integer :: na, nb
integer :: i, k, j, l, e
integer(bit_kind) :: deti(N_int), detk(N_int), detj(N_int), detl(N_int)
double precision :: phase
integer :: degree, h1, h2, p1, p2
integer :: exc(0:2,2)
integer :: list(N_int*bit_kind_size), n_elements
na = n_det_alpha_unique
nb = n_det_beta_unique
A_ik(:,:) = 0.d0
B_jl(:,:) = 0.d0
! -----------------------------------------------------------------------------------------------------------
! Diagonal part
! -----------------------------------------------------------------------------------------------------------
do i = 1, na
deti(1:N_int) = psi_det_alpha_unique(1:N_int,i)
list = 0
call bitstring_to_list(deti, list, n_elements, N_int)
do e = 1, elec_alpha_num
A_ik(i,i) += mo_one_e_integrals(list(e),list(e))
enddo
enddo
do j = 1, nb
detj(1:N_int) = psi_det_beta_unique(1:N_int,j)
list = 0
call bitstring_to_list(detj, list, n_elements, N_int)
do e = 1, elec_beta_num
B_jl(j,j) += mo_one_e_integrals(list(e),list(e))
enddo
enddo
! -----------------------------------------------------------------------------------------------------------
! -----------------------------------------------------------------------------------------------------------
! -----------------------------------------------------------------------------------------------------------
! degree = 1
! -----------------------------------------------------------------------------------------------------------
do i = 1, na
deti(1:N_int) = psi_det_alpha_unique(1:N_int,i)
list = 0
call bitstring_to_list(deti, list, n_elements, N_int)
do k = 1, na
detk(1:N_int) = psi_det_alpha_unique(1:N_int,k)
call get_excitation_degree_spin(deti, detk, degree, N_int)
if(degree .eq. 1) then
exc = 0
call get_single_excitation_spin(deti, detk, exc, phase, N_int)
call decode_exc_spin(exc, h1, p1, h2, p2)
A_ik(i,k) += phase * mo_one_e_integrals(h1,p1)
!A_ik(i,k) += phase * mo_one_e_integrals(p1,h1)
endif
enddo
enddo
do j = 1, nb
detj(1:N_int) = psi_det_beta_unique(1:N_int,j)
list = 0
call bitstring_to_list(detj, list, n_elements, N_int)
do l = 1, nb
detl(1:N_int) = psi_det_beta_unique(1:N_int,l)
call get_excitation_degree_spin(detj, detl, degree, N_int)
if(degree .eq. 1) then
exc = 0
call get_single_excitation_spin(detj, detl, exc, phase, N_int)
call decode_exc_spin(exc, h1, p1, h2, p2)
B_jl(j,l) += phase * mo_one_e_integrals(h1,p1)
!B_jl(j,l) += phase * mo_one_e_integrals(p1,h1)
endif
enddo
enddo
! -----------------------------------------------------------------------------------------------------------
! -----------------------------------------------------------------------------------------------------------
return
end subroutine const_AB
subroutine const_2b(T)
USE bitmasks
implicit none
double precision, external :: get_two_e_integral
double precision, intent(out) :: T(n_det_alpha_unique,n_det_alpha_unique,n_det_beta_unique,n_det_beta_unique )
integer :: na, nb
integer :: i, k, j, l
integer(bit_kind) :: psi_ij(N_int,2), psi_kl(N_int,2)
double precision :: phase
integer :: degree, h1, h2, p1, p2, s1, s2, e1, e2
integer :: ii, jj
integer :: exc(0:2,2,2)
integer :: occ(N_int*bit_kind_size,2), n_occ_alpha
double precision :: two_body_fact
na = n_det_alpha_unique
nb = n_det_beta_unique
T(:,:,:,:) = 0.d0
! -----------------------------------------------------------------------------------------------------------------
do i = 1, na
psi_ij(1:N_int,1) = psi_det_alpha_unique(1:N_int,i)
do j = 1, nb
psi_ij(1:N_int,2) = psi_det_beta_unique(1:N_int,j)
call bitstring_to_list(psi_ij(1,1), occ(1,1), n_occ_alpha, N_int)
call bitstring_to_list(psi_ij(1,2), occ(1,2), n_occ_alpha, N_int)
do k = 1, na
psi_kl(1:N_int,1) = psi_det_alpha_unique(1:N_int,k)
do l = 1, nb
psi_kl(1:N_int,2) = psi_det_beta_unique(1:N_int,l)
call get_excitation_degree(psi_ij, psi_kl, degree, N_int)
two_body_fact = 0.d0
if(degree .eq. 2) then
call get_double_excitation(psi_ij, psi_kl, exc, phase, N_int)
call decode_exc(exc, degree, h1, p1, h2, p2, s1, s2)
select case(s1+s2)
case(2,4)
two_body_fact += phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map)
two_body_fact -= phase * get_two_e_integral(h1, h2, p2, p1, mo_integrals_map)
case(3)
two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map)
two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map)
end select
else if(degree .eq. 1) then
call get_single_excitation(psi_ij, psi_kl, exc, phase, N_int)
call decode_exc(exc, degree, h1, p1, h2, p2, s1, s2)
select case(s1)
case(1)
do ii = 1, elec_alpha_num
p2 = occ(ii,1)
h2 = p2
two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map)
two_body_fact -= 0.5d0 * phase * get_two_e_integral(h1, h2, p2, p1, mo_integrals_map)
two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map)
two_body_fact -= 0.5d0 * phase * get_two_e_integral(h2, h1, p1, p2, mo_integrals_map)
enddo
do ii = 1, elec_beta_num
p2 = occ(ii,2)
h2 = p2
two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map)
two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map)
enddo
case(2)
do ii = 1, elec_alpha_num
p2 = occ(ii,1)
h2 = p2
two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map)
two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map)
enddo
do ii = 1, elec_beta_num
p2 = occ(ii,2)
h2 = p2
two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map)
two_body_fact -= 0.5d0 * phase * get_two_e_integral(h1, h2, p2, p1, mo_integrals_map)
two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map)
two_body_fact -= 0.5d0 * phase * get_two_e_integral(h2, h1, p1, p2, mo_integrals_map)
enddo
end select
else if(degree .eq. 0) then
do ii = 1, elec_alpha_num
e1 = occ(ii,1)
do jj = 1, elec_alpha_num
e2 = occ(jj,1)
two_body_fact += 0.5d0 * get_two_e_integral(e1, e2, e1, e2, mo_integrals_map)
two_body_fact -= 0.5d0 * get_two_e_integral(e1, e2, e2, e1, mo_integrals_map)
enddo
do jj = 1, elec_beta_num
e2 = occ(jj,2)
two_body_fact += 0.5d0 * get_two_e_integral(e1, e2, e1, e2, mo_integrals_map)
two_body_fact += 0.5d0 * get_two_e_integral(e2, e1, e2, e1, mo_integrals_map)
enddo
enddo
do ii = 1, elec_beta_num
e1 = occ(ii,2)
do jj = 1, elec_beta_num
e2 = occ(jj,2)
two_body_fact += 0.5d0 * get_two_e_integral(e1, e2, e1, e2, mo_integrals_map)
two_body_fact -= 0.5d0 * get_two_e_integral(e1, e2, e2, e1, mo_integrals_map)
enddo
enddo
end if
T(i,k,j,l) = two_body_fact
enddo
enddo
enddo
enddo
! -----------------------------------------------------------------------------------------------------------------
return
end subroutine const_2b
subroutine const_Hdiag(n_toselect, Uref, Vref, numalpha_toselect, numbeta_toselect &
, A_ik, B_jl, T, Hdiag)
implicit none
integer, intent(in) :: n_toselect
integer, intent(in) :: numalpha_toselect(n_toselect), numbeta_toselect(n_toselect)
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) :: A_ik(n_det_alpha_unique,n_det_alpha_unique)
double precision, intent(in) :: B_jl(n_det_beta_unique ,n_det_beta_unique )
double precision, intent(in) :: T(n_det_alpha_unique,n_det_alpha_unique,n_det_beta_unique,n_det_beta_unique )
double precision, intent(out) :: Hdiag(n_toselect)
integer :: na, nb
integer :: i, j, k, l, ii, jj, n
na = n_det_alpha_unique
nb = n_det_beta_unique
Hdiag(:) = 0.d0
do n = 1, n_toselect
ii = numalpha_toselect(n)
jj = numbeta_toselect (n)
do i = 1, na
do k = 1, na
Hdiag(n) += Uref(i,ii) * Uref(k,ii) * A_ik(i,k)
enddo
enddo
do j = 1, nb
do l = 1, nb
Hdiag(n) += Vref(j,jj) * Vref(l,jj) * B_jl(j,l)
enddo
enddo
do i = 1, na
do k = 1, na
do j = 1, nb
do l = 1, nb
Hdiag(n) += Uref(i,ii) * Uref(k,ii) * Vref(j,jj) * Vref(l,jj) * T(i,k,j,l)
enddo
enddo
enddo
enddo
enddo
return
end subroutine const_Hdiag