mirror of
https://gitlab.com/scemama/qp_plugins_scemama.git
synced 2024-11-07 22:53:42 +01:00
115 lines
3.3 KiB
Fortran
115 lines
3.3 KiB
Fortran
! !!!
|
|
subroutine RSVD( i_state, low_rank, PowerIt_max, nb_oversamp, URSVD, DRSVD, VtRSVD )
|
|
! !!!
|
|
BEGIN_DOC
|
|
! standard RSVD for a prefixed rank
|
|
END_DOC
|
|
! !!!
|
|
implicit none
|
|
include 'constants.include.F'
|
|
! !!!
|
|
integer, intent(in) :: i_state, low_rank, PowerIt_max, nb_oversamp
|
|
double precision, intent(out) :: URSVD(n_det_alpha_unique,low_rank), DRSVD(low_rank), VtRSVD(low_rank,n_det_beta_unique)
|
|
! !!!
|
|
integer :: i, j, k, l, PowerIt, m, n
|
|
double precision, allocatable :: r1(:,:), Q(:,:), P(:,:), B(:,:)
|
|
double precision, allocatable :: UB(:,:), D(:), Vt(:,:), U(:,:)
|
|
! !!!
|
|
m = n_det_alpha_unique
|
|
n = n_det_beta_unique
|
|
! !!!
|
|
! !!! !!! !!! !!! !!! !!! !!! !!! !!! !!! !
|
|
! !!!
|
|
allocate( Q(m, low_rank+nb_oversamp) )
|
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,k,l,r1)
|
|
allocate( r1(N_det,2) )
|
|
!$OMP DO
|
|
do l = 1, low_rank+nb_oversamp
|
|
Q(:,l) = 0.d0
|
|
call random_number(r1)
|
|
r1(:,1) = dsqrt(-2.d0*dlog(r1(:,1)))
|
|
r1(:,1) = r1(:,1) * dcos(dtwo_pi*r1(:,2))
|
|
do k = 1, N_det
|
|
i = psi_bilinear_matrix_rows(k)
|
|
Q(i,l) = Q(i,l) + psi_bilinear_matrix_values(k,i_state) * r1(k,1)
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
deallocate(r1)
|
|
!$OMP END PARALLEL
|
|
! !!!
|
|
call ortho_qr(Q, size(Q,1), m, low_rank+nb_oversamp)
|
|
! !!!
|
|
! power scheme
|
|
! !!!
|
|
allocate( P(n, low_rank+nb_oversamp) )
|
|
do PowerIt = 1, PowerIt_max
|
|
! !!!
|
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l)
|
|
!$OMP DO
|
|
do l = 1, low_rank+nb_oversamp
|
|
P(:,l) = 0.d0
|
|
do k = 1, N_det
|
|
i = psi_bilinear_matrix_rows(k)
|
|
j = psi_bilinear_matrix_columns(k)
|
|
P(j,l) = P(j,l) + psi_bilinear_matrix_values(k,i_state) * Q(i,l)
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP END PARALLEL
|
|
! !!!
|
|
call ortho_qr(P, size(P,1), n, low_rank+nb_oversamp)
|
|
! !!!
|
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l)
|
|
!$OMP DO
|
|
do l = 1, low_rank+nb_oversamp
|
|
Q(:,l) = 0.d0
|
|
do k = 1 , N_det
|
|
i = psi_bilinear_matrix_rows(k)
|
|
j = psi_bilinear_matrix_columns(k)
|
|
Q(i,l) = Q(i,l) + psi_bilinear_matrix_values(k,i_state) * P(j,l)
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP END PARALLEL
|
|
! !!!
|
|
call ortho_qr(Q, size(Q,1), m, low_rank+nb_oversamp)
|
|
! !!!
|
|
enddo
|
|
deallocate(P)
|
|
! !!!
|
|
allocate( B(low_rank+nb_oversamp,n) )
|
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l)
|
|
!$OMP DO
|
|
do l = 1, low_rank+nb_oversamp
|
|
B(l,:) = 0.d0
|
|
do k = 1, N_det
|
|
i = psi_bilinear_matrix_rows(k)
|
|
j = psi_bilinear_matrix_columns(k)
|
|
B(l,j) = B(l,j) + psi_bilinear_matrix_values(k,i_state) * Q(i,l)
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP END PARALLEL
|
|
! !!!
|
|
allocate( UB(low_rank+nb_oversamp,low_rank+nb_oversamp), D(low_rank+nb_oversamp), Vt(low_rank+nb_oversamp,n) )
|
|
call svd_s( B, size(B,1), UB, size(UB,1), D, Vt, size(Vt,1), low_rank+nb_oversamp, n)
|
|
deallocate(B)
|
|
allocate( U(m,low_rank+nb_oversamp) )
|
|
call dgemm('N', 'N', m, low_rank+nb_oversamp, low_rank+nb_oversamp, 1.d0, Q, size(Q,1), UB, size(UB,1), 0.d0, U, size(U,1))
|
|
deallocate( Q,UB )
|
|
! !!!
|
|
do l = 1, low_rank
|
|
DRSVD(l) = D(l)
|
|
do i = 1, m
|
|
URSVD(i,l) = U(i,l)
|
|
enddo
|
|
do i = 1, n
|
|
VtRSVD(l,i) = Vt(l,i)
|
|
enddo
|
|
enddo
|
|
! !!!
|
|
return
|
|
! !!!
|
|
end
|