10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-25 22:52:15 +02:00

dirty - corrected mrcepa/mrcc final PT2

This commit is contained in:
Yann Garniron 2016-06-02 12:47:35 +02:00
parent 42dc213725
commit acc8a8bb7e
4 changed files with 96 additions and 17 deletions

View File

@ -22,7 +22,7 @@ IRPF90_FLAGS : --ninja --align=32
# 0 : Deactivate
#
[OPTION]
MODE : DEBUG ; [ OPT | PROFILE | DEBUG ] : Chooses the section below
MODE : OPT ; [ OPT | PROFILE | DEBUG ] : Chooses the section below
CACHE : 1 ; Enable cache_compile.py
OPENMP : 1 ; Append OpenMP flags

View File

@ -31,5 +31,11 @@ s.set_perturbation("epstein_nesbet_2x2")
s.unset_openmp()
print s
s = H_apply_zmq("mrcepa_PT2")
s.energy = "psi_ref_energy_diagonalized"
s.set_perturbation("epstein_nesbet_2x2")
s.unset_openmp()
print s
END_SHELL

View File

@ -73,6 +73,7 @@ END_PROVIDER
BEGIN_PROVIDER [ double precision, lambda_mrcc, (N_states,psi_det_size) ]
&BEGIN_PROVIDER [ integer, lambda_mrcc_pt2, (0:psi_det_size) ]
&BEGIN_PROVIDER [ integer, lambda_mrcc_pt3, (0:psi_det_size) ]
implicit none
BEGIN_DOC
! cm/<Psi_0|H|D_m> or perturbative 1/Delta_E(m)
@ -81,7 +82,7 @@ END_PROVIDER
double precision :: ihpsi_current(N_states)
integer :: i_pert_count
double precision :: hii, lambda_pert
integer :: N_lambda_mrcc_pt2
integer :: N_lambda_mrcc_pt2, N_lambda_mrcc_pt3
integer :: histo(200), j
histo = 0
@ -92,7 +93,9 @@ END_PROVIDER
i_pert_count = 0
lambda_mrcc = 0.d0
N_lambda_mrcc_pt2 = 0
N_lambda_mrcc_pt3 = 0
lambda_mrcc_pt2(0) = 0
lambda_mrcc_pt3(0) = 0
do i=1,N_det_non_ref
call i_h_psi(psi_non_ref(1,1,i), psi_ref, psi_ref_coef, N_int, N_det_ref,&
@ -111,15 +114,21 @@ END_PROVIDER
N_lambda_mrcc_pt2 += 1
lambda_mrcc_pt2(N_lambda_mrcc_pt2) = i
endif
else
if (lambda_mrcc_pt3(N_lambda_mrcc_pt3) /= i) then
N_lambda_mrcc_pt3 += 1
lambda_mrcc_pt3(N_lambda_mrcc_pt3) = i
endif
endif
enddo
enddo
lambda_mrcc_pt2(0) = N_lambda_mrcc_pt2
lambda_mrcc_pt3(0) = N_lambda_mrcc_pt3
end if
print*,'N_det_non_ref = ',N_det_non_ref
print*,'Number of ignored determinants = ',i_pert_count
print*,'psi_coef_ref_ratio = ',psi_ref_coef(2,1)/psi_ref_coef(1,1)
print*,'lambda max = ',maxval(dabs(lambda_mrcc))
print*,'Number of ignored determinants = ',i_pert_count
END_PROVIDER

View File

@ -64,8 +64,10 @@ subroutine run_pt2(N_st,energy)
double precision, allocatable :: pt2(:), norm_pert(:), H_pert_diag(:)
integer, intent(in) :: N_st
double precision, intent(in) :: energy(N_st)
double precision :: pt3(N_st)
allocate (pt2(N_st), norm_pert(N_st),H_pert_diag(N_st))
pt2 = 0.d0
pt3 = 0d0
!if(lambda_mrcc_pt2(0) == 0) return
print*,'Last iteration only to compute the PT2'
@ -85,29 +87,91 @@ subroutine run_pt2(N_st,energy)
! enddo
! SOFT_TOUCH N_det_generators psi_det_generators psi_coef_generators ci_eigenvectors_dressed ci_eigenvectors_s2_dressed ci_electronic_energy_dressed
N_det_generators = lambda_mrcc_pt2(0) + N_det_cas
do i=1,N_det_cas
do k=1,N_int
psi_det_generators(k,1,i) = psi_ref(k,1,i)
psi_det_generators(k,2,i) = psi_ref(k,2,i)
enddo
do k=1,N_st
psi_coef_generators(i,k) = psi_ref_coef(i,k)
enddo
enddo
do i=N_det_cas+1,N_det_generators
j = lambda_mrcc_pt2(i - N_det_cas)
!
! N_det_generators = lambda_mrcc_pt2(0) + N_det_cas
! do i=1,N_det_cas
! do k=1,N_int
! psi_det_generators(k,1,i) = psi_ref(k,1,i)
! psi_det_generators(k,2,i) = psi_ref(k,2,i)
! enddo
! do k=1,N_st
! psi_coef_generators(i,k) = psi_ref_coef(i,k)
! enddo
! enddo
! do i=N_det_cas+1,N_det_generators
! j = lambda_mrcc_pt2(i - N_det_cas)
! do k=1,N_int
! psi_det_generators(k,1,i) = psi_non_ref(k,1,j)
! psi_det_generators(k,2,i) = psi_non_ref(k,2,j)
! enddo
! do k=1,N_st
! psi_coef_generators(i,k) = psi_non_ref_coef(j,k)
! enddo
! enddo
! SOFT_TOUCH N_det_generators psi_det_generators psi_coef_generators ci_eigenvectors_dressed ci_eigenvectors_s2_dressed ci_electronic_energy_dressed
N_det_generators = lambda_mrcc_pt3(0) + N_det_ref
N_det_selectors = lambda_mrcc_pt3(0) + N_det_ref
psi_det_generators(:,:,:N_det_ref) = psi_ref(:,:,:N_det_ref)
psi_selectors(:,:,:N_det_ref) = psi_ref(:,:,:N_det_ref)
psi_coef_generators(:N_det_ref,:) = psi_ref_coef(:N_det_ref,:)
psi_selectors_coef(:N_det_ref,:) = psi_ref_coef(:N_det_ref,:)
do i=N_det_ref+1,N_det_generators
j = lambda_mrcc_pt3(i-N_det_ref)
do k=1,N_int
psi_det_generators(k,1,i) = psi_non_ref(k,1,j)
psi_det_generators(k,2,i) = psi_non_ref(k,2,j)
psi_selectors(k,1,i) = psi_non_ref(k,1,j)
psi_selectors(k,2,i) = psi_non_ref(k,2,j)
enddo
do k=1,N_st
psi_coef_generators(i,k) = psi_non_ref_coef(j,k)
psi_selectors_coef(i,k) = psi_non_ref_coef(j,k)
enddo
enddo
SOFT_TOUCH N_det_generators psi_det_generators psi_coef_generators ci_eigenvectors_dressed ci_eigenvectors_s2_dressed ci_electronic_energy_dressed
SOFT_TOUCH N_det_selectors psi_selectors_coef psi_selectors N_det_generators psi_det_generators psi_coef_generators ci_eigenvectors_dressed ci_eigenvectors_s2_dressed ci_electronic_energy_dressed
SOFT_TOUCH psi_ref_coef_diagonalized psi_ref_energy_diagonalized! psi_coef_energy_diagonalized
call H_apply_mrcepa_PT2(pt2, norm_pert, H_pert_diag, N_st)
N_det_generators = N_det_non_ref + N_det_ref
N_det_selectors = N_det_non_ref + N_det_ref
psi_det_generators(:,:,:N_det_ref) = psi_ref(:,:,:N_det_ref)
psi_selectors(:,:,:N_det_ref) = psi_ref(:,:,:N_det_ref)
psi_coef_generators(:N_det_ref,:) = psi_ref_coef(:N_det_ref,:)
psi_selectors_coef(:N_det_ref,:) = psi_ref_coef(:N_det_ref,:)
do i=N_det_ref+1,N_det_generators
j = i-N_det_ref
do k=1,N_int
psi_det_generators(k,1,i) = psi_non_ref(k,1,j)
psi_det_generators(k,2,i) = psi_non_ref(k,2,j)
psi_selectors(k,1,i) = psi_non_ref(k,1,j)
psi_selectors(k,2,i) = psi_non_ref(k,2,j)
enddo
do k=1,N_st
psi_coef_generators(i,k) = psi_non_ref_coef(j,k)
psi_selectors_coef(i,k) = psi_non_ref_coef(j,k)
enddo
enddo
SOFT_TOUCH N_det_selectors psi_selectors_coef psi_selectors N_det_generators psi_det_generators psi_coef_generators ci_eigenvectors_dressed ci_eigenvectors_s2_dressed ci_electronic_energy_dressed
SOFT_TOUCH psi_ref_coef_diagonalized psi_ref_energy_diagonalized! psi_coef_energy_diagonalized
call H_apply_mrcepa_PT2(pt3, norm_pert, H_pert_diag, N_st)
call H_apply_mrcc_PT2(pt2, norm_pert, H_pert_diag, N_st)
!!!!!!!!!!!!!!!!
print *, "2-3 :",pt2, pt3
print *, lambda_mrcc_pt3(0), N_det, N_det_ref, psi_coef(1,1), psi_ref_coef(1,1)
pt2 = pt2 - pt3
print *, 'Final step'
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states