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_4bTSVD_v0.irp.f
2021-11-02 16:18:07 +01:00

306 lines
9.1 KiB
Fortran

program delta_4bTSVD_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, l
integer :: na, nb, nsvd
integer :: n_TSVD, La, lla, Lb, llb
integer :: n_tmp
double precision, allocatable :: U(:,:), V(:,:)
double precision, allocatable :: delta0_det(:), delta1_det(:)
double precision, allocatable :: delta_O(:,:), delta_D(:), delta_A(:,:), delta_B(:,:)
double precision, allocatable :: delta_tmp(:)
double precision :: t1, t2
na = n_det_alpha_unique
nb = n_det_beta_unique
call ezfio_get_spindeterminants_n_svd_coefs(nsvd)
allocate( u(na,nsvd) , v(nb,nsvd) )
print *, ' start read SVD vectors'
call wall_time(t1)
call ezfio_get_spindeterminants_psi_svd_alpha(U)
call ezfio_get_spindeterminants_psi_svd_beta (V)
call wall_time(t2)
print *, ' end read SVD vectors after (min)', (t2-t1)/60.d0
allocate( delta0_det(N_det) )
call ezfio_get_dmc_dress_dmc_delta_h(delta0_det)
! ---------------------------------------------
! ---------------------------------------------
! parameters
!
!n_TSVD = 100
read *, n_TSVD
!La = nsvd
!lla = 3
read *, La
read *, lla
!Lb = nsvd
!llb = 3
read *, Lb
read *, llb
print *, " --- delta_4bTSVD_v0 ---"
print *, " n_TSVD =", n_TSVD
print *, " La =", La
print *, " lla =", lla
print *, " Lb =", Lb
print *, " llb =", llb
print *, " tot =", (La-n_TSVD)*lla + &
(Lb-n_TSVD)*llb + &
n_TSVD * n_TSVD + &
nsvd - n_TSVD
!
! ---------------------------------------------
! ---------------------------------------------
! ---------------------------------------------
n_tmp = n_TSVD + 1
allocate( delta_O(1:n_TSVD,1:n_TSVD) )
allocate( delta_A(n_tmp:La,1:lla) )
allocate( delta_B(1:llb,n_tmp:Lb) )
allocate( delta_D(n_tmp:nsvd) )
call obtenir_deltaSVD_4b( nsvd, n_TSVD, La, lla, Lb, llb, U, V, delta0_det &
, delta_O, delta_A, delta_B, delta_D )
deallocate( delta0_det )
allocate( delta1_det(N_det) )
delta1_det(1:N_det) = 0.d0
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(nsvd, n_TSVD, n_tmp, La, lla, Lb, llb, N_det, U, V &
!$OMP , delta1_det, delta_O, delta_D, delta_A, delta_B &
!$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, n_TSVD
delta_tmp(ii) = delta_tmp(ii) + U(ii_i,l) * V(ii_j,l) * delta_D(l)
do k = 1, n_TSVD
delta_tmp(ii) = delta_tmp(ii) + U(ii_i,k) * V(ii_j,l) * delta_O(k,l)
enddo
enddo
do l = 1, lla
do k = n_tmp, La
delta_tmp(ii) = delta_tmp(ii) + U(ii_i,k) * V(ii_j,l) * delta_A(k,l)
enddo
enddo
do l = n_tmp, Lb
do k = 1, llb
delta_tmp(ii) = delta_tmp(ii) + U(ii_i,k) * V(ii_j,l) * delta_B(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
deallocate( U , V )
deallocate( delta_O, delta_A, delta_B, delta_D )
call ezfio_set_dmc_dress_dmc_delta_h(delta1_det)
deallocate( 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_4b( nsvd, n_TSVD, La, lla, Lb, llb, U, V, delta0_det &
, delta_O, delta_A, delta_B, delta_D )
implicit none
integer, intent(in) :: nsvd, n_TSVD, La, lla, Lb, llb
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) :: delta_O(1:n_TSVD,1:n_TSVD)
double precision, intent(out) :: delta_A(n_TSVD+1:La,1:lla)
double precision, intent(out) :: delta_B(1:llb,n_TSVD+1:Lb)
double precision, intent(out) :: delta_D(n_TSVD+1:nsvd)
integer :: n_tmp
integer :: k, l, ii, ii_i, ii_j
double precision, allocatable :: delta_tmp(:,:), delta1d_tmp(:)
n_tmp = n_TSVD + 1
! --------------------------------------------------------------------------------
!
delta_O(1:n_TSVD,1:n_TSVD) = 0.d0
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(n_TSVD, N_det, U, V, delta0_det, delta_O &
!$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) &
!$OMP PRIVATE(k, l, ii, ii_i, ii_j, delta_tmp)
allocate( delta_tmp(1:n_TSVD,1:n_TSVD) )
delta_tmp(1:n_TSVD,1:n_TSVD) = 0.d0
!$OMP DO
do l = 1, n_TSVD
do k = 1, n_TSVD
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, n_TSVD
do k = 1, n_TSVD
delta_O(k,l) = delta_O(k,l) + delta_tmp(k,l)
enddo
enddo
!$OMP END CRITICAL
deallocate(delta_tmp)
!$OMP END PARALLEL
!
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
!
delta_D(n_tmp:nsvd) = 0.d0
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(nsvd, n_tmp, N_det, U, V, delta0_det, delta_D &
!$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) &
!$OMP PRIVATE(l, ii, ii_i, ii_j, delta1d_tmp)
allocate( delta1d_tmp(n_tmp:nsvd) )
delta1d_tmp(n_tmp:nsvd) = 0.d0
!$OMP DO
do l = n_tmp, nsvd
do ii = 1, N_det
ii_i = psi_bilinear_matrix_rows (ii)
ii_j = psi_bilinear_matrix_columns(ii)
delta1d_tmp(l) = delta1d_tmp(l) + U(ii_i,l) * V(ii_j,l) * delta0_det(ii)
enddo
enddo
!$OMP END DO
!$OMP CRITICAL
do l = n_tmp, nsvd
delta_D(l) = delta_D(l) + delta1d_tmp(l)
enddo
!$OMP END CRITICAL
deallocate(delta1d_tmp)
!$OMP END PARALLEL
!
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
!
delta_A(n_tmp:La,1:lla) = 0.d0
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(n_tmp, La, lla, N_det, U, V, delta0_det, delta_A &
!$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) &
!$OMP PRIVATE(k, l, ii, ii_i, ii_j, delta_tmp)
allocate( delta_tmp(n_tmp:La,1:lla) )
delta_tmp(n_tmp:La,1:lla) = 0.d0
!$OMP DO
do l = 1, lla
do k = n_tmp, La
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, lla
do k = n_tmp, La
delta_A(k,l) = delta_A(k,l) + delta_tmp(k,l)
enddo
enddo
!$OMP END CRITICAL
deallocate(delta_tmp)
!$OMP END PARALLEL
!
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
!
delta_B(1:llb,n_tmp:Lb) = 0.d0
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(n_tmp, Lb, llb, N_det, U, V, delta0_det, delta_B &
!$OMP , psi_bilinear_matrix_rows, psi_bilinear_matrix_columns) &
!$OMP PRIVATE(k, l, ii, ii_i, ii_j, delta_tmp)
allocate( delta_tmp(1:llb,n_tmp:Lb) )
delta_tmp(1:llb,n_tmp:Lb) = 0.d0
!$OMP DO
do l = n_tmp, Lb
do k = 1, llb
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 = n_tmp, Lb
do k = 1, llb
delta_B(k,l) = delta_B(k,l) + delta_tmp(k,l)
enddo
enddo
!$OMP END CRITICAL
deallocate(delta_tmp)
!$OMP END PARALLEL
!
! --------------------------------------------------------------------------------
return
end subroutine obtenir_deltaSVD_4b
! _______________________________________________________________________________________________________
! _______________________________________________________________________________________________________