1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-06-02 03:15:25 +02:00
qp_plugins_scemama/devel/svdwf/print_ij_H_kl_v0.irp.f
2021-11-02 16:18:07 +01:00

226 lines
6.3 KiB
Fortran

program print_ij_H_kl_v0
implicit none
BEGIN_DOC
! perturbative approach to build psi_postsvd
! without OMP
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
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(:)
integer :: cantor_pairing_ij, cantor_pairing_new
integer, allocatable :: cantor_pairing(:), cantor_pairing_tmp(:)
double precision :: t_beg, t_end
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 *, ' 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 cpu_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 cpu_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 )
! get < u_i v_j | H | u_k v_l > for all vectors
call const_ij_H_kl(Uref, Vref)
deallocate( Uref, Vref )
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 *, " Execution avec 1 threads"
print *, " total elapsed time (min) = ", W_tot_time/60.d0
print *, ' ___________________________________________________________________'
end
subroutine const_ij_H_kl(Uref, Vref)
implicit none
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)
integer(bit_kind) :: det1(N_int,2)
integer(bit_kind) :: det2(N_int,2)
integer :: degree, na, nb
integer :: i, j, k, l, ii, jj, kk, ll
double precision :: h12, h_abgd, sum_i, sum_j, sum_k, sum_l
double precision :: t1, t2
print *, ""
print *, " start const_ij_H_kl"
call wall_time(t1)
na = n_det_alpha_unique
nb = n_det_beta_unique
do ii = 1, na
do jj = 1, nb
do kk = 1, na
do ll = 1, nb
sum_l = 0.d0
do l = 1, nb
det2(:,2) = psi_det_beta_unique(:,l)
sum_j = 0.d0
do j = 1, nb
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
sum_k = 0.d0
do k = 1, na
det2(:,1) = psi_det_alpha_unique(:,k)
sum_i = 0.d0
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
sum_i = sum_i + h12 * Uref(i,ii)
enddo
sum_k = sum_k + sum_i * Uref(k,kk)
enddo
sum_j = sum_j + sum_k * Vref(j,jj)
enddo
sum_l = sum_l + sum_j * Vref(l,ll)
enddo
h_abgd = sum_l
write(7070, '(4(I4,2X),(F15.8))') ii, jj, kk, ll, h_abgd
enddo
enddo
enddo
enddo
call wall_time(t2)
print *, " end const_ij_H_kl after (min) ", (t2-t1)/60.
print *, ""
return
end subroutine const_ij_H_kl