10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-02 03:15:29 +02:00
quantum_package/plugins/mrcepa0/mrcepa0_general.irp.f

254 lines
7.3 KiB
Fortran
Raw Normal View History

2016-04-07 17:54:13 +02:00
2016-03-04 16:52:46 +01:00
2016-04-07 17:54:13 +02:00
subroutine run(N_st,energy)
2016-03-04 16:52:46 +01:00
implicit none
2016-04-07 17:54:13 +02:00
integer, intent(in) :: N_st
double precision, intent(out) :: energy(N_st)
2016-04-11 17:42:15 +02:00
integer :: i,j
2016-03-04 16:52:46 +01:00
double precision :: E_new, E_old, delta_e
2016-04-07 17:54:13 +02:00
integer :: iteration
2016-11-10 14:42:41 +01:00
double precision :: E_past(4)
2016-04-07 17:54:13 +02:00
integer :: n_it_mrcc_max
double precision :: thresh_mrcc
2016-11-10 14:42:41 +01:00
double precision, allocatable :: lambda(:)
allocate (lambda(N_states))
2016-05-20 09:44:22 +02:00
2016-10-07 19:21:04 +02:00
thresh_mrcc = thresh_dressed_ci
n_it_mrcc_max = n_it_max_dressed_ci
2016-04-07 17:54:13 +02:00
2016-07-06 15:43:21 +02:00
if(n_it_mrcc_max == 1) then
2016-11-21 23:31:28 +01:00
do j=1,N_states
2016-04-11 17:42:15 +02:00
do i=1,N_det
psi_coef(i,j) = CI_eigenvectors_dressed(i,j)
enddo
enddo
SOFT_TOUCH psi_coef ci_energy_dressed
call write_double(6,ci_energy_dressed(1),"Final MRCC energy")
2016-10-14 12:40:29 +02:00
call ezfio_set_mrcepa0_energy(ci_energy_dressed(1))
2016-04-11 17:42:15 +02:00
call save_wavefunction
else
E_new = 0.d0
delta_E = 1.d0
iteration = 0
lambda = 1.d0
do while (delta_E > thresh_mrcc)
iteration += 1
2016-11-21 21:25:38 +01:00
print *, '==============================================='
2017-09-14 17:20:42 +02:00
print *, 'Iteration', iteration, '/', n_it_mrcc_max
2016-11-21 21:25:38 +01:00
print *, '==============================================='
2016-04-11 17:42:15 +02:00
print *, ''
E_old = mrcc_e0_denominator(1) !sum(ci_energy_dressed(1:N_states))
2016-11-20 22:55:10 +01:00
do i=1,N_st
2017-09-14 17:20:42 +02:00
call write_double(6,ci_energy_dressed(i),"Energy")
2016-11-20 22:55:10 +01:00
enddo
2016-04-11 17:42:15 +02:00
call diagonalize_ci_dressed(lambda)
E_new = mrcc_e0_denominator(1) !sum(ci_energy_dressed(1:N_states))
2017-12-06 16:16:44 +01:00
2017-06-02 01:15:29 +02:00
! if (.true.) then
! provide delta_ij_mrcc_pouet
! endif
2016-11-20 22:55:10 +01:00
delta_E = (E_new - E_old)/dble(N_states)
2016-11-21 21:25:38 +01:00
print *, ''
call write_double(6,thresh_mrcc,"thresh_mrcc")
2016-11-20 22:55:10 +01:00
call write_double(6,delta_E,"delta_E")
delta_E = dabs(delta_E)
2016-04-11 17:42:15 +02:00
call save_wavefunction
2016-10-07 19:21:04 +02:00
call ezfio_set_mrcepa0_energy(ci_energy_dressed(1))
if (iteration >= n_it_mrcc_max) then
2016-04-11 17:42:15 +02:00
exit
endif
enddo
2017-09-14 17:20:42 +02:00
call write_double(6,ci_energy_dressed(1),"Final energy")
2016-04-11 17:42:15 +02:00
endif
2016-11-14 19:43:07 +01:00
energy(1:N_st) = ci_energy_dressed(1:N_st)
2016-03-04 16:52:46 +01:00
end
2016-09-17 23:33:06 +02:00
subroutine print_cas_coefs
implicit none
integer :: i,j
print *, 'CAS'
print *, '==='
do i=1,N_det_cas
2016-11-09 15:42:27 +01:00
print *, (psi_cas_coef(i,j), j=1,N_states)
2016-09-17 23:33:06 +02:00
call debug_det(psi_cas(1,1,i),N_int)
enddo
call write_double(6,ci_energy(1),"Initial CI energy")
end
subroutine run_pt2_old(N_st,energy)
2016-04-07 17:54:13 +02:00
implicit none
integer :: i,j,k
integer, intent(in) :: N_st
double precision, intent(in) :: energy(N_st)
2016-09-17 23:33:06 +02:00
double precision :: pt2_redundant(N_st), pt2(N_st)
double precision :: norm_pert(N_st),H_pert_diag(N_st)
pt2_redundant = 0.d0
pt2 = 0d0
2016-04-07 17:54:13 +02:00
!if(lambda_mrcc_pt2(0) == 0) return
print*,'Last iteration only to compute the PT2'
2016-09-17 23:33:06 +02:00
print * ,'Computing the redundant PT2 contribution'
if (mrmode == 1) then
N_det_generators = lambda_mrcc_kept(0)
N_det_selectors = lambda_mrcc_kept(0)
do i=1,N_det_generators
j = lambda_mrcc_kept(i)
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
2016-05-04 18:23:05 +02:00
enddo
2016-09-17 23:33:06 +02:00
else
N_det_generators = N_det_non_ref
N_det_selectors = N_det_non_ref
do i=1,N_det_generators
j = i
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
2016-05-04 18:23:05 +02:00
enddo
2016-09-17 23:33:06 +02:00
endif
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
2016-09-17 23:33:06 +02:00
SOFT_TOUCH psi_ref_coef_diagonalized psi_ref_energy_diagonalized
call H_apply_mrcepa_PT2(pt2_redundant, norm_pert, H_pert_diag, N_st)
2016-09-17 23:33:06 +02:00
print * ,'Computing the remaining contribution'
2016-11-10 14:42:41 +01:00
threshold_selectors = max(threshold_selectors,threshold_selectors_pt2)
threshold_generators = max(threshold_generators,threshold_generators_pt2)
2016-09-17 23:33:06 +02:00
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
2016-05-04 18:23:05 +02:00
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)
2016-05-04 18:23:05 +02:00
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)
2016-05-04 18:23:05 +02:00
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
2016-09-17 23:33:06 +02:00
SOFT_TOUCH psi_ref_coef_diagonalized psi_ref_energy_diagonalized
2016-09-17 23:33:06 +02:00
call H_apply_mrcepa_PT2(pt2, norm_pert, H_pert_diag, N_st)
2016-09-17 23:33:06 +02:00
print *, "Redundant PT2 :",pt2_redundant
print *, "Full PT2 :",pt2
print *, lambda_mrcc_kept(0), N_det, N_det_ref, psi_coef(1,1), psi_ref_coef(1,1)
pt2 = pt2 - pt2_redundant
2016-04-07 17:54:13 +02:00
print *, 'Final step'
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states
print *, 'PT2 = ', pt2
print *, 'E = ', energy
print *, 'E+PT2 = ', energy+pt2
print *, '-----'
2016-10-07 19:21:04 +02:00
call ezfio_set_mrcepa0_energy_pt2(energy(1)+pt2(1))
2016-04-07 17:54:13 +02:00
end
2016-09-17 23:33:06 +02:00
subroutine run_pt2(N_st,energy)
implicit none
integer :: i,j,k
integer, intent(in) :: N_st
double precision, intent(in) :: energy(N_st)
double precision :: pt2(N_st)
double precision :: norm_pert(N_st),H_pert_diag(N_st)
pt2 = 0d0
!if(lambda_mrcc_pt2(0) == 0) return
print*,'Last iteration only to compute the PT2'
N_det_generators = N_det_cas
N_det_selectors = N_det_non_ref
2016-04-07 17:54:13 +02:00
2016-09-17 23:33:06 +02:00
do i=1,N_det_generators
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=1,N_det
do k=1,N_int
psi_selectors(k,1,i) = psi_det_sorted(k,1,i)
psi_selectors(k,2,i) = psi_det_sorted(k,2,i)
enddo
do k=1,N_st
psi_selectors_coef(i,k) = psi_coef_sorted(i,k)
enddo
2016-03-04 16:52:46 +01:00
enddo
2016-04-07 17:54:13 +02:00
2016-09-17 23:33:06 +02:00
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
2016-03-04 16:52:46 +01:00
2016-09-17 23:33:06 +02:00
call H_apply_mrcepa_PT2(pt2, norm_pert, H_pert_diag, N_st)
! call ezfio_set_full_ci_energy_pt2(energy+pt2)
print *, 'Final step'
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states
print *, 'PT2 = ', pt2
print *, 'E = ', energy
print *, 'E+PT2 = ', energy+pt2
print *, '-----'
2016-10-07 19:21:04 +02:00
call ezfio_set_mrcepa0_energy_pt2(energy(1)+pt2(1))
2016-09-17 23:33:06 +02:00
end
2016-03-04 16:52:46 +01:00