10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-13 17:43:55 +01:00
quantum_package/plugins/Full_CI_ZMQ/selection.irp.f

107 lines
3.0 KiB
Fortran
Raw Normal View History

2016-09-05 17:16:09 +02:00
use bitmasks
2016-09-08 15:22:09 +02:00
double precision function integral8(i,j,k,l)
2016-09-05 17:16:09 +02:00
implicit none
2016-09-08 15:22:09 +02:00
integer, intent(in) :: i,j,k,l
double precision, external :: get_mo_bielec_integral
2016-09-05 17:16:09 +02:00
2016-09-08 15:22:09 +02:00
integral8 = get_mo_bielec_integral(i,j,k,l,mo_integrals_map)
end function
2016-09-05 17:16:09 +02:00
2016-09-08 17:34:56 +02:00
BEGIN_PROVIDER [ integer(1), psi_phasemask, (N_int*bit_kind_size, 2, N_det)]
2016-09-05 17:16:09 +02:00
use bitmasks
implicit none
integer :: i
do i=1, N_det
call get_mask_phase(psi_det_sorted(1,1,i), psi_phasemask(1,1,i))
end do
END_PROVIDER
subroutine assert(cond, msg)
character(*), intent(in) :: msg
logical, intent(in) :: cond
if(.not. cond) then
print *, "assert fail: "//msg
stop
end if
end subroutine
subroutine get_mask_phase(det, phasemask)
use bitmasks
implicit none
integer(bit_kind), intent(in) :: det(N_int, 2)
2016-09-08 17:34:56 +02:00
integer(1), intent(out) :: phasemask(N_int*bit_kind_size, 2)
2016-09-05 17:16:09 +02:00
integer :: s, ni, i
logical :: change
2016-09-08 17:34:56 +02:00
phasemask = 0_1
2016-09-05 17:16:09 +02:00
do s=1,2
change = .false.
do ni=1,N_int
do i=0,bit_kind_size-1
if(BTEST(det(ni, s), i)) change = .not. change
2016-09-08 17:34:56 +02:00
if(change) phasemask((ni-1)*bit_kind_size + i + 1, s) = 1_1
2016-09-05 17:16:09 +02:00
end do
end do
end do
end subroutine
subroutine select_connected(i_generator,E0,pt2,b)
use bitmasks
use selection_types
implicit none
integer, intent(in) :: i_generator
type(selection_buffer), intent(inout) :: b
double precision, intent(inout) :: pt2(N_states)
integer :: k,l
double precision, intent(in) :: E0(N_states)
integer(bit_kind) :: hole_mask(N_int,2), particle_mask(N_int,2)
double precision :: fock_diag_tmp(2,mo_tot_num+1)
call build_fock_tmp(fock_diag_tmp,psi_det_generators(1,1,i_generator),N_int)
do l=1,N_generators_bitmask
do k=1,N_int
hole_mask(k,1) = iand(generators_bitmask(k,1,s_hole,l), psi_det_generators(k,1,i_generator))
hole_mask(k,2) = iand(generators_bitmask(k,2,s_hole,l), psi_det_generators(k,2,i_generator))
particle_mask(k,1) = iand(generators_bitmask(k,1,s_part,l), not(psi_det_generators(k,1,i_generator)) )
particle_mask(k,2) = iand(generators_bitmask(k,2,s_part,l), not(psi_det_generators(k,2,i_generator)) )
enddo
2016-09-05 17:18:01 +02:00
call select_doubles(i_generator,hole_mask,particle_mask,fock_diag_tmp,E0,pt2,b)
2016-09-05 17:16:09 +02:00
call select_singles(i_generator,hole_mask,particle_mask,fock_diag_tmp,E0,pt2,b)
enddo
end subroutine
double precision function get_phase_bi(phasemask, s1, s2, h1, p1, h2, p2)
use bitmasks
implicit none
2016-09-08 17:34:56 +02:00
integer(1), intent(in) :: phasemask(N_int*bit_kind_size, 2)
2016-09-05 17:16:09 +02:00
integer, intent(in) :: s1, s2, h1, h2, p1, p2
logical :: change
2016-09-08 17:34:56 +02:00
integer(1) :: np
2016-09-08 10:12:28 +02:00
double precision, parameter :: res(0:1) = (/1d0, -1d0/)
2016-09-08 17:34:56 +02:00
np = phasemask(h1,s1) + phasemask(p1,s1) + phasemask(h2,s2) + phasemask(p2,s2)
if(p1 < h1) np = np + 1_1
if(p2 < h2) np = np + 1_1
if(s1 == s2 .and. max(h1, p1) > min(h2, p2)) np = np + 1_1
get_phase_bi = res(iand(np,1_1))
2016-09-05 17:16:09 +02:00
end subroutine