program buildpsi_diagSVDit_Anthony_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, 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(:,:), H(:,:,:,:) double precision, allocatable :: psi_postsvd(:,:), 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) 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) i_state = 1 ! --------------------------------------------------------------------------------------- ! 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 ) ! --------------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------------- ! numerote vectors ! Truncated rank n_selected = 20 print*, ' initial psi space rank = ', n_selected ! check SVD error err0 = 0.d0 do j = 1, n_det_beta_unique 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 ) ! --------------------------------------------------------------------------------------- !________________________________________________________________________________________________________ ! ! loop over SVD iterations !________________________________________________________________________________________________________ E0_old = 0.d0 tol_energy = 1.d0 it_svd = 0 it_svd_max = 100 allocate( H(n_det_beta_unique,n_det_beta_unique,n_det_beta_unique,n_det_beta_unique) ) allocate( psi_postsvd(n_det_beta_unique,n_det_beta_unique) ) 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 double precision :: norm norm = 0.d0 do j = 1, n_selected norm = norm + Dref(j)*Dref(j) enddo Dref = Dref / dsqrt(norm) print *, '' print *, '' print *, '' print *, '-- Compute H --' call const_H_uv(Uref, Vref, H) ! H0(i,j) = < u_i v_j | H | u_i v_j > ! E0 = < psi_0 | H | psi_0 > E0 = 0.d0 do j = 1, n_selected do i = 1, n_selected E0 = E0 + Dref(i) * H(i,i,j,j) * Dref(j) enddo enddo E0_av = E0 + nuclear_repulsion print *,' E0 (avant SVD) =', E0_av double precision, allocatable :: eigval0(:) double precision, allocatable :: eigvec0(:,:,:) double precision, allocatable :: H_tmp(:,:,:,:) allocate( H_tmp(n_selected,n_selected,n_selected,n_selected) ) do l=1,n_selected do k=1,n_selected do j=1,n_selected do i=1,n_selected H_tmp(i,j,k,l) = H(i,j,k,l) enddo enddo enddo enddo allocate( eigval0(n_selected**2),eigvec0(n_selected,n_selected,n_selected**2)) eigvec0 = 0.d0 print *, ' --- Diag post-SVD --- ' call lapack_diag(eigval0, eigvec0, H_tmp, n_selected**2, n_selected**2) E0_postsvd = eigval0(1)+nuclear_repulsion print*, ' postsvd energy = ', E0_postsvd deallocate(H_tmp, eigval0) print *, ' --- SVD --- ' Dref = 0.d0 call perform_newpostSVD(n_selected, eigvec0(1,1,1), Uref, Vref, Dref) deallocate(eigvec0) print *, ' --- Compute H --- ' call const_H_uv(Uref, Vref, H) ! H0(i,j) = < u_i v_j | H | u_i v_j > ! E0 = < psi_0 | H | psi_0 > E0 = 0.d0 norm = 0.d0 do j = 1, n_det_beta_unique do i = 1, n_det_beta_unique E0 = E0 + Dref(i) * H(i,i,j,j) * Dref(j) enddo norm = norm + Dref(j)*Dref(j) enddo E0_ap = E0 + nuclear_repulsion print *,' E0 (apres SVD) =', E0_ap !print *,' norm =', norm print *, ' --- Perturbation --- ' psi_postsvd = 0.d0 do i=1,n_selected psi_postsvd(i,i) = Dref(i) enddo !do j=1,n_selected ! do i=n_selected+1,n_det_beta_unique ! print *, i,j, H(i,j,i,j) ! enddo !enddo !do j=n_selected+1,n_det_beta_unique ! do i=1,n_selected ! print *, i,j, H(i,j,i,j) ! enddo !enddo !do j=n_selected+1,n_det_beta_unique ! do i=n_selected+1,n_det_beta_unique ! print *, i,j, H(i,j,i,j) ! enddo !enddo Ept2 = 0.d0 do j=1,n_selected do i=n_selected+1,n_det_beta_unique ctmp = 0.d0 do l=1,n_selected do k=1,n_selected ctmp = ctmp + H(k,l,i,j) * psi_postsvd(k,l) enddo enddo psi_postsvd(i,j) = ctmp / (E0_ap - (H(i,j,i,j)+nuclear_repulsion) ) Ept2 += ctmp*ctmp / (E0_ap - (H(i,j,i,j)+nuclear_repulsion) ) enddo enddo do j=n_selected+1,n_det_beta_unique do i=1,n_selected ctmp = 0.d0 do l=1,n_selected do k=1,n_selected ctmp = ctmp + H(k,l,i,j) * psi_postsvd(k,l) enddo enddo psi_postsvd(i,j) = ctmp / (E0_ap - (H(i,j,i,j)+nuclear_repulsion) ) Ept2 += ctmp*ctmp / (E0_ap - (H(i,j,i,j)+nuclear_repulsion) ) enddo enddo do j=n_selected+1,n_det_beta_unique do i=n_selected+1,n_det_beta_unique ctmp = 0.d0 do l=1,n_selected do k=1,n_selected ctmp = ctmp + H(k,l,i,j) * psi_postsvd(k,l) enddo enddo psi_postsvd(i,j) = ctmp / (E0_ap - (H(i,j,i,j)+nuclear_repulsion) ) Ept2 += ctmp*ctmp / (E0_ap - (H(i,j,i,j)+nuclear_repulsion) ) enddo enddo E0pt2 = E0_ap + Ept2 print *, ' perturb energy = ', E0pt2, Ept2 tol_energy = 100.d0 * dabs(E0pt2-E0_old) / dabs(E0pt2) E0_old = E0pt2 print *, ' --- SVD --- ' call perform_newpostSVD(n_det_beta_unique, psi_postsvd, Uref, Vref, Dref) write(22,'(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( H, psi_postsvd ) 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 perform_newpostSVD(n_selected, psi_postsvd, Uref, Vref, Dref) USE OMP_LIB integer, intent(in) :: n_selected double precision, intent(in) :: psi_postsvd(n_selected,n_selected) 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 j = 1, n_selected do i = 1, n_selected S_mat(i,j) = psi_postsvd(i,j) 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 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, 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 subroutine const_H_uv(Uref, Vref, H) USE OMP_LIB implicit none 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) :: H(n_det_beta_unique,n_det_beta_unique, n_det_beta_unique,n_det_beta_unique) 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 integer :: p,q,r,s double precision :: h12, x double precision, allocatable :: H0(:,:,:,:) double precision, allocatable :: H1(:,:,:,:) allocate( H0(n_det_alpha_unique,n_det_beta_unique, n_det_alpha_unique, n_det_beta_unique) ) allocate( H1(n_det_alpha_unique,n_det_beta_unique, n_det_alpha_unique, n_det_beta_unique) ) H0(:,:,:,:) = 0.d0 !$OMP PARALLEL DEFAULT(NONE) & !$OMP PRIVATE(p,q,r,s,i,j,k,l,det1,det2,degree) & !$OMP SHARED(n_det_alpha_unique,n_det_beta_unique,psi_det_alpha_unique,psi_det_beta_unique, & !$OMP N_int,Uref,Vref,H0,H1,H) !$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(det2, det1, N_int, H0(k,l,i,j) ) ! !!! enddo enddo enddo enddo !$OMP END DO !$OMP SINGLE H1 = 0.d0 !$OMP END SINGLE !$OMP DO do s = 1, n_det_beta_unique do l = 1, n_det_beta_unique do k = 1, n_det_alpha_unique do j = 1, n_det_beta_unique do i = 1, n_det_alpha_unique H1(i,j,k,s) = H1(i,j,k,s) + H0(i,j,k,l) * Vref(l,s) enddo enddo enddo enddo enddo !$OMP END DO !$OMP SINGLE H0 = 0.d0 !$OMP END SINGLE !$OMP DO do s = 1, n_det_beta_unique do r = 1, n_det_beta_unique do k = 1, n_det_alpha_unique do j = 1, n_det_beta_unique do i = 1, n_det_alpha_unique H0(i,j,r,s) = H0(i,j,r,s) + H1(i,j,k,s) * Uref(k,r) enddo enddo enddo enddo enddo !$OMP END DO !$OMP SINGLE H1 = 0.d0 !$OMP END SINGLE !$OMP DO do s = 1, n_det_beta_unique do j = 1, n_det_beta_unique do r = 1, n_det_alpha_unique do q = 1, n_det_beta_unique do i = 1, n_det_alpha_unique H1(i,q,r,s) = H1(i,q,r,s) + H0(i,j,r,s) * Vref(j,q) enddo enddo enddo enddo enddo !$OMP END DO !$OMP SINGLE H = 0.d0 !$OMP END SINGLE !$OMP DO do s = 1, n_det_beta_unique do r = 1, n_det_beta_unique do q = 1, n_det_beta_unique do p = 1, n_det_beta_unique do i = 1, n_det_alpha_unique H(p,q,r,s) = H(p,q,r,s) + H1(i,q,r,s) * Uref(i,p) enddo enddo enddo enddo enddo !$OMP END DO !$OMP END PARALLEL deallocate(H1,H0) return end