10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-05 13:43:57 +01:00
quantum_package/plugins/mrcc_sto/mrcc_sto.irp.f

241 lines
7.5 KiB
Fortran
Raw Normal View History

2018-02-19 17:15:59 +01:00
program mrcc_sto
implicit none
BEGIN_DOC
! TODO
END_DOC
call dress_zmq()
2018-02-27 18:43:07 +01:00
call ezfio_set_mrcc_sto_energy(ci_energy_dressed(1))
end
2018-02-20 15:51:53 +01:00
2018-02-20 17:34:51 +01:00
BEGIN_PROVIDER [ double precision, hij_cache_, (N_det,Nproc) ]
&BEGIN_PROVIDER [ double precision, sij_cache_, (N_det,Nproc) ]
&BEGIN_PROVIDER [ double precision, dIa_hla_, (N_states,N_det,Nproc) ]
&BEGIN_PROVIDER [ double precision, dIa_sla_, (N_states,N_det,Nproc) ]
2018-02-26 11:33:32 +01:00
&BEGIN_PROVIDER [ integer, excs_ , (0:2,2,2,N_det,Nproc) ]
&BEGIN_PROVIDER [ double precision, phases_, (N_det, Nproc) ]
BEGIN_DOC
2018-03-05 17:04:26 +01:00
! temporay arrays for dress_with_alpha_buffer. Avoids reallocation.
2018-02-20 17:34:51 +01:00
END_DOC
END_PROVIDER
2018-02-22 15:29:20 +01:00
2018-03-02 15:29:58 +01:00
subroutine dress_with_alpha_buffer(delta_ij_loc, i_gen, minilist, det_minilist, n_minilist, alpha, iproc)
2018-02-19 17:15:59 +01:00
use bitmasks
implicit none
2018-02-20 15:51:53 +01:00
BEGIN_DOC
!delta_ij_loc(:,:,1) : dressing column for H
!delta_ij_loc(:,:,2) : dressing column for S2
!minilist : indices of determinants connected to alpha ( in psi_det_sorted )
!n_minilist : size of minilist
!alpha : alpha determinant
END_DOC
2018-02-22 15:29:20 +01:00
integer(bit_kind), intent(in) :: alpha(N_int,2), det_minilist(N_int, 2, n_minilist)
2018-03-02 15:29:58 +01:00
integer,intent(in) :: minilist(n_minilist), n_minilist, iproc, i_gen
2018-02-19 17:15:59 +01:00
double precision, intent(inout) :: delta_ij_loc(N_states,N_det,2)
integer :: i,j,k,l,m
integer :: degree1, degree2, degree
double precision :: hIk, hla, hIl, sla, dIk(N_states), dka(N_states), dIa(N_states), hka
double precision :: phase, phase2
integer :: exc(0:2,2,2)
integer :: h1,h2,p1,p2,s1,s2
2018-02-20 15:51:53 +01:00
integer(bit_kind) :: tmp_det(N_int,2), ctrl
2018-02-19 17:15:59 +01:00
integer :: i_state, k_sd, l_sd, m_sd, ll_sd, i_I
double precision :: Delta_E_inv(N_states)
double precision :: sdress, hdress
2018-02-20 15:51:53 +01:00
logical :: ok, ok2
2018-02-26 11:33:32 +01:00
integer :: canbediamond
2018-02-20 17:34:51 +01:00
PROVIDE mo_class
2018-02-19 17:15:59 +01:00
2018-02-20 17:34:51 +01:00
2018-02-20 15:51:53 +01:00
if(n_minilist == 1) return
2018-02-22 13:41:11 +01:00
do i=1,n_minilist
2018-02-23 14:02:51 +01:00
if(idx_non_ref_rev(minilist(i)) == 0) return
2018-02-22 13:41:11 +01:00
end do
2018-02-20 15:51:53 +01:00
2018-02-19 17:15:59 +01:00
if (perturbative_triples) then
PROVIDE one_anhil fock_virt_total fock_core_inactive_total one_creat
endif
2018-03-01 11:35:00 +01:00
2018-02-26 11:33:32 +01:00
canbediamond = 0
2018-02-19 17:15:59 +01:00
do l_sd=1,n_minilist
2018-03-01 11:35:00 +01:00
call get_excitation(det_minilist(1,1,l_sd),alpha,exc,degree1,phase,N_int)
call decode_exc(exc,degree1,h1,p1,h2,p2,s1,s2)
2018-02-26 11:33:32 +01:00
2018-03-01 11:35:00 +01:00
ok = (mo_class(h1)(1:1) == 'A' .or. mo_class(h1)(1:1) == 'I') .and. &
(mo_class(p1)(1:1) == 'A' .or. mo_class(p1)(1:1) == 'V')
if(ok .and. degree1 == 2) then
ok = (mo_class(h2)(1:1) == 'A' .or. mo_class(h2)(1:1) == 'I') .and. &
(mo_class(p2)(1:1) == 'A' .or. mo_class(p2)(1:1) == 'V')
end if
if(ok) then
canbediamond += 1
excs_(:,:,:,l_sd,iproc) = exc(:,:,:)
phases_(l_sd, iproc) = phase
else
phases_(l_sd, iproc) = 0d0
end if
!call i_h_j(alpha,det_minilist(1,1,l_sd),N_int,hij_cache_(l_sd,iproc))
!call get_s2(alpha,det_minilist(1,1,l_sd),N_int,sij_cache_(l_sd,iproc))
call i_h_j_s2(alpha,det_minilist(1,1,l_sd),N_int,hij_cache_(l_sd,iproc), sij_cache_(l_sd,iproc))
2018-02-19 17:15:59 +01:00
enddo
2018-02-26 11:33:32 +01:00
if(canbediamond <= 1) return
2018-02-19 17:15:59 +01:00
do i_I=1,N_det_ref
call get_excitation_degree(alpha,psi_ref(1,1,i_I),degree1,N_int)
if (degree1 > 4) then
cycle
endif
do i_state=1,N_states
dIa(i_state) = 0.d0
enddo
2018-02-23 14:02:51 +01:00
do k_sd=1,n_minilist
2018-02-26 11:33:32 +01:00
if(phases_(k_sd,iproc) == 0d0) cycle
2018-02-23 14:02:51 +01:00
call get_excitation_degree(psi_ref(1,1,i_I),det_minilist(1,1,k_sd),degree,N_int)
2018-02-19 17:15:59 +01:00
if (degree > 2) then
cycle
endif
2018-02-26 11:33:32 +01:00
!call get_excitation(det_minilist(1,1,k_sd),alpha,exc,degree2,phase,N_int)
phase = phases_(k_sd, iproc)
exc(:,:,:) = excs_(:,:,:,k_sd,iproc)
degree2 = exc(0,1,1) + exc(0,1,2)
2018-02-19 17:15:59 +01:00
call apply_excitation(psi_ref(1,1,i_I), exc, tmp_det, ok, N_int)
2018-02-23 14:02:51 +01:00
if((.not. ok) .and. (.not. perturbative_triples)) cycle
2018-02-19 17:15:59 +01:00
do i_state=1,N_states
dka(i_state) = 0.d0
enddo
2018-02-23 14:02:51 +01:00
ok2 = .false.
!do i_state=1,N_states
! !if(dka(i_state) == 0) cycle
! dIk(i_state) = dij(i_I, idx_non_ref_rev(minilist(k_sd)), i_state)
! if(dIk(i_state) /= 0d0) then
! ok2 = .true.
! endif
!enddo
!if(.not. ok2) cycle
2018-02-19 17:15:59 +01:00
if (ok) then
2018-02-23 14:02:51 +01:00
phase2 = 0d0
do l_sd=k_sd+1,n_minilist
2018-02-26 11:33:32 +01:00
if(phases_(l_sd, iproc) == 0d0) cycle
2018-02-23 14:02:51 +01:00
call get_excitation_degree(tmp_det,det_minilist(1,1,l_sd),degree,N_int)
2018-02-19 17:15:59 +01:00
if (degree == 0) then
do i_state=1,N_states
2018-02-23 14:02:51 +01:00
dIk(i_state) = dij(i_I, idx_non_ref_rev(minilist(k_sd)), i_state)
if(dIk(i_state) /= 0d0) then
if(phase2 == 0d0) call get_excitation(psi_ref(1,1,i_I),det_minilist(1,1,l_sd),exc,degree,phase2,N_int)
dka(i_state) = dij(i_I, idx_non_ref_rev(minilist(l_sd)), i_state) * phase * phase2
end if
end do
!call get_excitation(psi_ref(1,1,i_I),det_minilist(1,1,l_sd),exc,degree,phase2,N_int)
!do i_state=1,N_states
! if(dIk(i_state) /= 0d0) dka(i_state) = dij(i_I, idx_non_ref_rev(minilist(l_sd)), i_state) * phase * phase2
!enddo
2018-02-19 17:15:59 +01:00
exit
2018-02-23 14:02:51 +01:00
2018-02-19 17:15:59 +01:00
endif
enddo
else if (perturbative_triples) then
2018-02-23 14:02:51 +01:00
hka = hij_cache_(k_sd,iproc)
if (dabs(hka) > 1.d-12) then
call get_delta_e_dyall_general_mp(psi_ref(1,1,i_I),alpha,Delta_E_inv)
do i_state=1,N_states
ASSERT (Delta_E_inv(i_state) < 0.d0)
dka(i_state) = hka / Delta_E_inv(i_state)
enddo
endif
2018-02-19 17:15:59 +01:00
endif
2018-02-23 14:02:51 +01:00
2018-02-19 17:15:59 +01:00
if (perturbative_triples.and. (degree2 == 1) ) then
call i_h_j(psi_ref(1,1,i_I),tmp_det,N_int,hka)
2018-02-23 14:02:51 +01:00
hka = hij_cache_(k_sd,iproc) - hka
2018-02-19 17:15:59 +01:00
if (dabs(hka) > 1.d-12) then
call get_delta_e_dyall_general_mp(psi_ref(1,1,i_I),alpha,Delta_E_inv)
do i_state=1,N_states
ASSERT (Delta_E_inv(i_state) < 0.d0)
dka(i_state) = hka / Delta_E_inv(i_state)
enddo
endif
endif
do i_state=1,N_states
dIa(i_state) = dIa(i_state) + dIk(i_state) * dka(i_state)
enddo
enddo
2018-02-20 15:51:53 +01:00
ok2 = .false.
do i_state=1,N_states
if(dIa(i_state) /= 0d0) ok2 = .true.
enddo
if(.not. ok2) cycle
2018-02-23 14:02:51 +01:00
do l_sd=1,n_minilist
k_sd = minilist(l_sd)
hla = hij_cache_(l_sd,iproc)
sla = sij_cache_(l_sd,iproc)
2018-02-19 17:15:59 +01:00
do i_state=1,N_states
2018-02-20 17:34:51 +01:00
hdress = dIa(i_state) * hla * psi_ref_coef(i_I,i_state)
sdress = dIa(i_state) * sla * psi_ref_coef(i_I,i_state)
2018-03-02 15:29:58 +01:00
!!!$OMP ATOMIC
2018-02-23 14:02:51 +01:00
delta_ij_loc(i_state,k_sd,1) += hdress
2018-03-02 15:29:58 +01:00
!!!$OMP ATOMIC
2018-02-23 14:02:51 +01:00
delta_ij_loc(i_state,k_sd,2) += sdress
2018-02-19 17:15:59 +01:00
enddo
enddo
enddo
end subroutine
!! TESTS MINILIST
subroutine test_minilist(minilist, n_minilist, alpha)
use bitmasks
implicit none
integer, intent(in) :: n_minilist
integer(bit_kind),intent(in) :: alpha(N_int, 2)
integer, intent(in) :: minilist(n_minilist)
integer :: a, i, deg
integer :: refc(N_det), testc(N_det)
refc = 0
testc = 0
do i=1,N_det
2018-02-26 11:33:32 +01:00
call get_excitation_degree(psi_det(1,1,i), alpha, deg, N_int)
if(deg <= 2) refc(i) = refc(i) + 1
end do
do i=1,n_minilist
2018-02-26 11:33:32 +01:00
call get_excitation_degree(psi_det(1,1,minilist(i)), alpha, deg, N_int)
if(deg <= 2) then
testc(minilist(i)) += 1
else
stop "NON LINKED IN MINILIST"
end if
end do
do i=1,N_det
if(refc(i) /= testc(i)) then
print *, "MINILIST FAIL ", sum(refc), sum(testc), n_minilist
exit
end if
end do
end subroutine