1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-06-02 03:15:25 +02:00
qp_plugins_scemama/devel/svdwf/delta_svd_v1.irp.f
2021-11-02 16:18:07 +01:00

222 lines
6.8 KiB
Fortran

program delta_svd_v1
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 :: deltaDet_ex(:), deltaDet_ap(:), deltaSVD_kl(:,:)
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( deltaSVD_kk(nsvd) )
!call obtenir_deltaSVD_kk(nsvd, U, V, deltaSVD_kk)
allocate( deltaSVD_kl(nsvd,nsvd) )
call obtenir_deltaSVD_kl(nsvd, U, V, deltaSVD_kl)
allocate( deltaDet_ex(N_det) , deltaDet_ap(N_det) )
do ii = 1, N_det
ii_i = psi_bilinear_matrix_rows (ii)
ii_j = psi_bilinear_matrix_columns(ii)
deltaDet_ex(ii) = 0.d0
deltaDet_ap(ii) = 0.d0
do k = 1, nsvd
deltaDet_ap(ii) = deltaDet_ap(ii) + U(ii_i,k) * V(ii_j,k) * deltaSVD_kl(k,k)
do l = 1, nsvd
deltaDet_ex(ii) = deltaDet_ex(ii) + U(ii_i,k) * V(ii_j,l) * deltaSVD_kl(k,l)
enddo
enddo
print *, deltaDet_ex(ii) , deltaDet_ap(ii)
enddo
!call ezfio_set_dmc_dress_dmc_delta_h(delta_Det_ex)
deallocate( U , S , V )
deallocate( deltaSVD_kl )
deallocate( deltaDet_ex , deltaDet_ap )
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, 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(out) :: deltaSVD_kk(nsvd)
integer :: i_state, na, nb
integer :: k, ii, jj, ii_i, ii_j, jj_i, jj_j
double precision :: tmp_k, tmp_ij, coef_jj
integer(bit_kind) :: det_ii(N_int,2), det_jj(N_int,2)
double precision :: hmono, heff, hderiv, hthree, htilde_ij, hij, delta_mat
PROVIDE scalar_mu_r_pot_physicist_mo deriv_mu_r_pot_physicist_mo
PROVIDE three_body_3_index three_body_3_index_exch_12 three_body_3_index_exch_13 three_body_3_index_exch_23
PROVIDE three_body_5_index three_body_5_index_exch_13 three_body_5_index_exch_32
PROVIDE three_body_4_index three_body_4_index_exch_12 three_body_4_index_exch_12_part
i_state = 1
na = n_det_alpha_unique
nb = n_det_beta_unique
deltaSVD_kk(:) = 0.d0
do k = 1, nsvd
tmp_k = 0.d0
do ii = 1, N_det
ii_i = psi_bilinear_matrix_rows (ii)
ii_j = psi_bilinear_matrix_columns(ii)
det_ii(:,1) = psi_det_alpha_unique(:,ii_i)
det_ii(:,2) = psi_det_beta_unique (:,ii_j)
!tmp_ij = < dup_i ddn_j | H_tilde - H | psi >
tmp_ij = 0.d0
do jj = 1, N_det
jj_i = psi_bilinear_matrix_rows (jj)
jj_j = psi_bilinear_matrix_columns(jj)
coef_jj = psi_bilinear_matrix_values(jj,i_state)
det_jj(:,1) = psi_det_alpha_unique(:,jj_i)
det_jj(:,2) = psi_det_beta_unique (:,jj_j)
call htilde_mat(det_ii, det_jj, hmono, heff, hderiv, hthree, htilde_ij)
call i_H_j(det_ii, det_jj, N_int, hij)
delta_mat = htilde_ij - hij
tmp_ij = tmp_ij + coef_jj * delta_mat
enddo
tmp_k = tmp_k + U(ii_i,k) * V(ii_j,k) * tmp_ij
enddo
deltaSVD_kk(k) = tmp_k
enddo
return
end subroutine obtenir_deltaSVD_kk
! _______________________________________________________________________________________________________
! _______________________________________________________________________________________________________
! _______________________________________________________________________________________________________
!
! [ 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, 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(out) :: deltaSVD_kl(nsvd,nsvd)
integer :: i_state, na, nb
integer :: k, l, ii, jj, ii_i, ii_j, jj_i, jj_j
double precision :: tmp_kl, tmp_ij, coef_jj
integer(bit_kind) :: det_ii(N_int,2), det_jj(N_int,2)
double precision :: hmono, heff, hderiv, hthree, htilde_ij, hij, delta_mat
PROVIDE scalar_mu_r_pot_physicist_mo deriv_mu_r_pot_physicist_mo
PROVIDE three_body_3_index three_body_3_index_exch_12 three_body_3_index_exch_13 three_body_3_index_exch_23
PROVIDE three_body_5_index three_body_5_index_exch_13 three_body_5_index_exch_32
PROVIDE three_body_4_index three_body_4_index_exch_12 three_body_4_index_exch_12_part
i_state = 1
na = n_det_alpha_unique
nb = n_det_beta_unique
! ! !
det_ii(:,1) = psi_det_alpha_unique(:,1)
det_ii(:,2) = psi_det_beta_unique (:,1)
det_jj(:,1) = psi_det_alpha_unique(:,1)
det_jj(:,2) = psi_det_beta_unique (:,1)
call htilde_mat(det_ii, det_jj, hmono, heff, hderiv, hthree, htilde_ij)
call i_H_j(det_ii, det_jj, N_int, hij)
! ! !
deltaSVD_kl(:,:) = 0.d0
do k = 1, nsvd
do l = 1, nsvd
tmp_kl = 0.d0
do ii = 1, N_det
ii_i = psi_bilinear_matrix_rows (ii)
ii_j = psi_bilinear_matrix_columns(ii)
det_ii(:,1) = psi_det_alpha_unique(:,ii_i)
det_ii(:,2) = psi_det_beta_unique (:,ii_j)
!tmp_ij = < dup_i ddn_j | H_tilde - H | psi >
tmp_ij = 0.d0
do jj = 1, N_det
jj_i = psi_bilinear_matrix_rows (jj)
jj_j = psi_bilinear_matrix_columns(jj)
coef_jj = psi_bilinear_matrix_values(jj,i_state)
det_jj(:,1) = psi_det_alpha_unique(:,jj_i)
det_jj(:,2) = psi_det_beta_unique (:,jj_j)
call htilde_mat(det_ii, det_jj, hmono, heff, hderiv, hthree, htilde_ij)
call i_H_j(det_ii, det_jj, N_int, hij)
delta_mat = htilde_ij - hij
tmp_ij = tmp_ij + coef_jj * delta_mat
enddo
tmp_kl = tmp_kl + U(ii_i,k) * V(ii_j,l) * tmp_ij
enddo
deltaSVD_kl(k,l) = tmp_kl
enddo
enddo
return
end subroutine obtenir_deltaSVD_kl
! _______________________________________________________________________________________________________
! _______________________________________________________________________________________________________