mirror of
https://gitlab.com/scemama/qp_plugins_scemama.git
synced 2024-12-23 04:43:38 +01:00
148 lines
3.9 KiB
Fortran
148 lines
3.9 KiB
Fortran
program delta_diagSVD_v0
|
|
|
|
BEGIN_DOC
|
|
! !!!
|
|
END_DOC
|
|
|
|
implicit none
|
|
|
|
read_wf = .True.
|
|
touch read_wf
|
|
|
|
call run()
|
|
|
|
end
|
|
|
|
|
|
|
|
subroutine run()
|
|
|
|
implicit none
|
|
|
|
integer :: ii, ii_i, ii_j, k
|
|
integer :: na, nb, nsvd, n_TSVD
|
|
double precision, allocatable :: U(:,:), V(:,:), S(:)
|
|
|
|
double precision, allocatable :: delta0_det(:), delta1_det(:), deltaSVD_kk(:)
|
|
double precision, allocatable :: delta_tmp(:)
|
|
|
|
na = n_det_alpha_unique
|
|
nb = n_det_beta_unique
|
|
|
|
! read SVD vectors
|
|
call ezfio_get_spindeterminants_n_svd_coefs(nsvd)
|
|
|
|
allocate( U(na,nsvd) , V(nb,nsvd) )
|
|
call ezfio_get_spindeterminants_psi_svd_alpha(U)
|
|
call ezfio_get_spindeterminants_psi_svd_beta (V)
|
|
|
|
allocate( delta0_det(N_det) )
|
|
call ezfio_get_dmc_dress_dmc_delta_h(delta0_det)
|
|
|
|
allocate( deltaSVD_kk(nsvd) )
|
|
deltaSVD_kk(1:nsvd) = 0.d0
|
|
call obtenir_deltaSVD_kk(nsvd, U, V, delta0_det, deltaSVD_kk)
|
|
deallocate( delta0_det )
|
|
|
|
|
|
! ---------------------------------------------
|
|
! parameters
|
|
|
|
n_TSVD = nsvd
|
|
|
|
print *, " --- delta_diagSVD_v0 ---"
|
|
print *, " n_TSVD =", n_TSVD
|
|
print *, " tot =", n_TSVD
|
|
|
|
!
|
|
! ---------------------------------------------
|
|
|
|
|
|
! [ delta ]_mn = \sum_{k} U_{m,k} delta_SVD_{k,k} V_{n,k}
|
|
allocate( delta1_det(N_det) )
|
|
delta1_det(:) = 0.d0
|
|
|
|
!$OMP PARALLEL DEFAULT(NONE) &
|
|
!$OMP SHARED(n_TSVD, N_det, U, V, deltaSVD_kk, delta1_det &
|
|
!$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) &
|
|
!$OMP PRIVATE(k, ii, ii_i, ii_j, delta_tmp)
|
|
allocate( delta_tmp(N_det) )
|
|
delta_tmp(:) = 0.d0
|
|
!$OMP DO
|
|
do ii = 1, N_det
|
|
ii_i = psi_bilinear_matrix_rows (ii)
|
|
ii_j = psi_bilinear_matrix_columns(ii)
|
|
do k = 1, n_TSVD
|
|
delta_tmp(ii) = delta_tmp(ii) + U(ii_i,k) * V(ii_j,k) * deltaSVD_kk(k)
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP CRITICAL
|
|
do ii = 1, N_det
|
|
delta1_det(ii) = delta1_det(ii) + delta_tmp(ii)
|
|
enddo
|
|
!$OMP END CRITICAL
|
|
deallocate( delta_tmp )
|
|
!$OMP END PARALLEL
|
|
|
|
deallocate( U , V )
|
|
deallocate( deltaSVD_kk )
|
|
|
|
call ezfio_set_dmc_dress_dmc_delta_h(delta1_det)
|
|
deallocate( delta1_det )
|
|
|
|
|
|
end
|
|
|
|
|
|
|
|
|
|
|
|
! _______________________________________________________________________________________________________
|
|
!
|
|
! [ delta_SVD ]_kk = < dup_SVD ddn_SVD | H_tilde - H | psi >
|
|
! = \sum_{i,j} U_{i,k} V_{j,k} < dup_i ddn_j | H_tilde - H | psi >
|
|
!
|
|
subroutine obtenir_deltaSVD_kk(nsvd, U, V, delta0_det, deltaSVD_kk)
|
|
|
|
implicit none
|
|
|
|
integer, intent(in) :: nsvd
|
|
double precision, intent(in) :: U(n_det_alpha_unique,nsvd), V(n_det_beta_unique,nsvd)
|
|
double precision, intent(in) :: delta0_det(N_det)
|
|
|
|
double precision, intent(out) :: deltaSVD_kk(nsvd)
|
|
|
|
integer :: k, ii, ii_i, ii_j
|
|
double precision, allocatable :: delta_tmp(:)
|
|
|
|
deltaSVD_kk(1:nsvd) = 0.d0
|
|
|
|
!$OMP PARALLEL DEFAULT(NONE) &
|
|
!$OMP SHARED(nsvd, N_det, U, V, delta0_det, deltaSVD_kk &
|
|
!$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) &
|
|
!$OMP PRIVATE(k, ii, ii_i, ii_j, delta_tmp)
|
|
allocate( delta_tmp(nsvd) )
|
|
delta_tmp(1:nsvd) = 0.d0
|
|
!$OMP DO
|
|
do k = 1, nsvd
|
|
do ii = 1, N_det
|
|
ii_i = psi_bilinear_matrix_rows (ii)
|
|
ii_j = psi_bilinear_matrix_columns(ii)
|
|
delta_tmp(k) = delta_tmp(k) + U(ii_i,k) * V(ii_j,k) * delta0_det(ii)
|
|
enddo
|
|
enddo
|
|
!$OMP END DO
|
|
!$OMP CRITICAL
|
|
do k = 1, nsvd
|
|
deltaSVD_kk(k) = deltaSVD_kk(k) + delta_tmp(k)
|
|
enddo
|
|
!$OMP END CRITICAL
|
|
deallocate(delta_tmp)
|
|
!$OMP END PARALLEL
|
|
|
|
return
|
|
end subroutine obtenir_deltaSVD_kk
|
|
! _______________________________________________________________________________________________________
|
|
! _______________________________________________________________________________________________________
|