mirror of
https://gitlab.com/scemama/qp_plugins_scemama.git
synced 2024-11-07 22:53:42 +01:00
98 lines
3.0 KiB
Fortran
98 lines
3.0 KiB
Fortran
|
|
program FSVD_trunc
|
|
implicit none
|
|
BEGIN_DOC
|
|
! study precision variation with truncated SVD
|
|
END_DOC
|
|
read_wf = .True.
|
|
TOUCH read_wf
|
|
! !!!
|
|
call run()
|
|
! !!!
|
|
end
|
|
|
|
|
|
subroutine run
|
|
|
|
implicit none
|
|
|
|
integer :: mm, nn, i_state, low_rank, lrank_min, lrank_max
|
|
integer :: i, j, k, l
|
|
double precision :: norm_psi, err_verif, err_tmp
|
|
double precision, allocatable :: U_FSVD(:,:), D_FSVD(:), Vt_FSVD(:,:), A_FSVD(:,:)
|
|
double precision, allocatable :: Uezfio(:,:,:), Dezfio(:,:), Vezfio(:,:,:)
|
|
|
|
i_state = 1
|
|
mm = n_det_alpha_unique
|
|
nn = n_det_beta_unique
|
|
|
|
print *, ' matrix dimensions:', mm,'x',nn
|
|
print *, ' N det:', N_det
|
|
|
|
allocate( A_FSVD(mm,nn), U_FSVD(mm,mm), D_FSVD(min(mm,nn)), Vt_FSVD(nn,nn) )
|
|
|
|
norm_psi = 0.d0
|
|
A_FSVD(:,:) = 0.d0
|
|
do k = 1, N_det
|
|
i = psi_bilinear_matrix_rows(k)
|
|
j = psi_bilinear_matrix_columns(k)
|
|
A_FSVD(i,j) = psi_bilinear_matrix_values(k,i_state)
|
|
norm_psi += psi_bilinear_matrix_values(k,i_state) * psi_bilinear_matrix_values(k,i_state)
|
|
enddo
|
|
|
|
call svd_s( A_FSVD, size(A_FSVD,1), U_FSVD, size(U_FSVD,1), D_FSVD, Vt_FSVD, size(Vt_FSVD,1), mm, nn)
|
|
print *, ' --- Full SVD: ok --- '
|
|
|
|
!lrank_min = 100
|
|
!lrank_max = nn
|
|
!do low_rank = lrank_min, lrank_max, 1
|
|
! err_verif = 0.d0
|
|
! do j = 1, nn
|
|
! do i = 1, mm
|
|
! err_tmp = 0.d0
|
|
! do l = 1, low_rank
|
|
! err_tmp = err_tmp + D_FSVD(l) * U_FSVD(i,l) * Vt_FSVD(l,j)
|
|
! enddo
|
|
! err_verif = err_verif + (A_FSVD(i,j)-err_tmp) * (A_FSVD(i,j)-err_tmp)
|
|
! enddo
|
|
! enddo
|
|
! print*, ' low rank = ', low_rank
|
|
! print*, ' err verif (%) = ', 100.d0 * dsqrt(err_verif/norm_psi)
|
|
!enddo
|
|
|
|
! ------------------------------------------------------------------------------------------------
|
|
! set to EZFIO for a fixed low rank
|
|
|
|
low_rank = min(mm,nn)
|
|
allocate( Uezfio(mm,low_rank,1), Dezfio(low_rank,1), Vezfio(nn,low_rank,1))
|
|
|
|
do l = 1, low_rank
|
|
Dezfio(l,1) = D_FSVD(l)
|
|
do j = 1, mm
|
|
Uezfio(j,l,1) = U_FSVD(j,l)
|
|
enddo
|
|
do j = 1, nn
|
|
Vezfio(j,l,1) = Vt_FSVD(l,j)
|
|
enddo
|
|
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(low_rank)
|
|
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 )
|
|
deallocate( U_FSVD, D_FSVD, Vt_FSVD )
|
|
|
|
end
|