program delta_FSVD_v2 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, l integer :: na, nb, nsvd double precision, allocatable :: U(:,:), V(:,:), S(:) double precision, allocatable :: delta0_det(:), delta1_det(:), deltaSVD_kl(:,:) 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), S(nsvd), V(nb,nsvd) ) call ezfio_get_spindeterminants_psi_svd_alpha(U) call ezfio_get_spindeterminants_psi_svd_beta (V) call ezfio_get_spindeterminants_psi_svd_coefs(S) allocate( delta0_det(N_det) ) call ezfio_get_dmc_dress_dmc_delta_h(delta0_det) allocate( deltaSVD_kl(nsvd,nsvd) ) deltaSVD_kl(:,:) = 0.d0 call obtenir_deltaSVD_kl(nsvd, U, V, delta0_det, deltaSVD_kl) ! [ delta ]_mn = \sum_{k,l} U_{m,k} delta_SVD_{k,l} V_{n,l} allocate( delta1_det(N_det) ) delta1_det(1:N_det) = 0.d0 !$OMP PARALLEL DEFAULT(NONE) & !$OMP SHARED(nsvd, N_det, U, V, deltaSVD_kl, delta1_det & !$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) & !$OMP PRIVATE(k, l, ii, ii_i, ii_j, delta_tmp) allocate( delta_tmp(N_det) ) delta_tmp(1:N_det) = 0.d0 !$OMP DO do ii = 1, N_det ii_i = psi_bilinear_matrix_rows (ii) ii_j = psi_bilinear_matrix_columns(ii) do l = 1, nsvd do k = 1, nsvd delta_tmp(ii) = delta_tmp(ii) + U(ii_i,k) * V(ii_j,l) * deltaSVD_kl(k,l) enddo 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 call ezfio_set_dmc_dress_dmc_delta_h(delta1_det) deallocate( U , S , V ) deallocate( deltaSVD_kl ) deallocate( delta0_det , delta1_det ) end ! _______________________________________________________________________________________________________ ! ! [ delta_SVD ]_kl = < dup_SVD ddn_SVD | H_tilde - H | psi > ! = \sum_{i,j} U_{i,k} V_{j,l} < dup_i ddn_j | H_tilde - H | psi > ! subroutine obtenir_deltaSVD_kl(nsvd, U, V, delta0_det, deltaSVD_kl) 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_kl(nsvd,nsvd) integer :: k, l, ii, ii_i, ii_j double precision, allocatable :: delta_tmp(:,:) deltaSVD_kl(1:nsvd,1:nsvd) = 0.d0 !$OMP PARALLEL DEFAULT(NONE) & !$OMP SHARED(nsvd, N_det, U, V, delta0_det, deltaSVD_kl & !$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) & !$OMP PRIVATE(k, l, ii, ii_i, ii_j, delta_tmp) allocate( delta_tmp(nsvd,nsvd) ) delta_tmp(1:nsvd,1:nsvd) = 0.d0 !$OMP DO do l = 1, nsvd 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,l) = delta_tmp(k,l) + U(ii_i,k) * V(ii_j,l) * delta0_det(ii) enddo enddo enddo !$OMP END DO !$OMP CRITICAL do l = 1, nsvd do k = 1, nsvd deltaSVD_kl(k,l) = deltaSVD_kl(k,l) + delta_tmp(k,l) enddo enddo !$OMP END CRITICAL deallocate(delta_tmp) !$OMP END PARALLEL return end subroutine obtenir_deltaSVD_kl ! _______________________________________________________________________________________________________ ! _______________________________________________________________________________________________________