10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-03 12:43:52 +01:00
quantum_package/plugins/Full_CI_ZMQ/selection_buffer.irp.f

78 lines
1.8 KiB
Fortran
Raw Normal View History

2016-09-01 14:43:13 +02:00
subroutine create_selection_buffer(N, siz, res)
use selection_types
implicit none
integer, intent(in) :: N, siz
type(selection_buffer), intent(out) :: res
allocate(res%det(N_int, 2, siz), res%val(siz))
res%val = 0d0
res%det = 0_8
res%N = N
res%mini = 0d0
res%cur = 0
end subroutine
subroutine add_to_selection_buffer(b, det, val)
use selection_types
implicit none
type(selection_buffer), intent(inout) :: b
integer(bit_kind), intent(in) :: det(N_int, 2)
double precision, intent(in) :: val
integer :: i
if(dabs(val) >= b%mini) then
b%cur += 1
2017-03-03 12:02:21 +01:00
b%det(1:N_int,1:2,b%cur) = det(1:N_int,1:2)
2016-09-01 14:43:13 +02:00
b%val(b%cur) = val
if(b%cur == size(b%val)) then
call sort_selection_buffer(b)
end if
end if
end subroutine
subroutine sort_selection_buffer(b)
use selection_types
implicit none
type(selection_buffer), intent(inout) :: b
2017-03-30 01:13:27 +02:00
double precision, allocatable:: absval(:)
2016-09-01 14:43:13 +02:00
integer, allocatable :: iorder(:)
2017-03-30 01:13:27 +02:00
double precision, pointer :: vals(:)
integer(bit_kind), pointer :: detmp(:,:,:)
2016-09-01 14:43:13 +02:00
integer :: i, nmwen
logical, external :: detEq
nmwen = min(b%N, b%cur)
2017-03-30 01:13:27 +02:00
allocate(iorder(b%cur), detmp(N_int, 2, size(b%det,3)), absval(b%cur), vals(size(b%val)))
2016-09-01 14:43:13 +02:00
absval = -dabs(b%val(:b%cur))
do i=1,b%cur
iorder(i) = i
end do
2017-03-30 01:13:27 +02:00
! Optimal for almost sorted data
2017-04-14 18:11:02 +02:00
! call sorted_dnumber(absval, b%cur, i)
! if (b%cur/i >
! call insertion_dsort(absval, iorder, b%cur)
call dsort(absval, iorder, b%cur)
2016-09-01 14:43:13 +02:00
do i=1, nmwen
2017-01-16 16:31:49 +01:00
detmp(1:N_int,1,i) = b%det(1:N_int,1,iorder(i))
detmp(1:N_int,2,i) = b%det(1:N_int,2,iorder(i))
2016-09-01 14:43:13 +02:00
vals(i) = b%val(iorder(i))
end do
2017-04-12 18:26:57 +02:00
do i=nmwen+1, size(vals)
vals(i) = 0.d0
enddo
2017-03-30 01:13:27 +02:00
deallocate(b%det, b%val)
b%det => detmp
b%val => vals
2016-09-01 14:43:13 +02:00
b%mini = max(b%mini,dabs(b%val(b%N)))
b%cur = nmwen
end subroutine