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

140 lines
3.3 KiB
Fortran
Raw Normal View History

2016-09-01 14:43:13 +02:00
2017-11-26 11:09:03 +01:00
subroutine create_selection_buffer(N, siz_, res)
2016-09-01 14:43:13 +02:00
use selection_types
implicit none
2017-11-26 11:09:03 +01:00
integer, intent(in) :: N, siz_
2016-09-01 14:43:13 +02:00
type(selection_buffer), intent(out) :: res
2017-11-26 11:09:03 +01:00
integer :: siz
siz = max(siz_,1)
2016-09-01 14:43:13 +02:00
allocate(res%det(N_int, 2, siz), res%val(siz))
2017-11-23 10:35:13 +01:00
res%val(:) = 0d0
res%det(:,:,:) = 0_8
2016-09-01 14:43:13 +02:00
res%N = N
res%mini = 0d0
res%cur = 0
end subroutine
2017-05-05 11:32:17 +02:00
subroutine delete_selection_buffer(b)
use selection_types
implicit none
type(selection_buffer), intent(inout) :: b
if (associated(b%det)) then
deallocate(b%det)
endif
if (associated(b%val)) then
deallocate(b%val)
endif
end
2016-09-01 14:43:13 +02:00
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(b%N > 0 .and. val <= b%mini) then
2016-09-01 14:43:13 +02:00
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
2017-05-05 11:32:17 +02:00
subroutine merge_selection_buffers(b1, b2)
use selection_types
implicit none
BEGIN_DOC
! Merges the selection buffers b1 and b2 into b2
END_DOC
2017-05-10 20:42:14 +02:00
type(selection_buffer), intent(inout) :: b1
2017-05-05 11:32:17 +02:00
type(selection_buffer), intent(inout) :: b2
integer(bit_kind), pointer :: detmp(:,:,:)
double precision, pointer :: val(:)
integer :: i, i1, i2, k, nmwen
2017-05-10 22:09:53 +02:00
if (b1%cur == 0) return
2017-05-10 20:42:14 +02:00
do while (b1%val(b1%cur) > b2%mini)
b1%cur = b1%cur-1
2017-05-10 21:33:53 +02:00
if (b1%cur == 0) then
return
endif
2017-05-10 20:42:14 +02:00
enddo
2017-05-10 21:50:25 +02:00
nmwen = min(b1%N, b1%cur+b2%cur)
allocate( val(size(b1%val)), detmp(N_int, 2, size(b1%det,3)) )
i1=1
i2=1
2017-05-05 11:32:17 +02:00
do i=1,nmwen
if ( (i1 > b1%cur).and.(i2 > b2%cur) ) then
exit
else if (i1 > b1%cur) then
val(i) = b2%val(i2)
detmp(1:N_int,1,i) = b2%det(1:N_int,1,i2)
detmp(1:N_int,2,i) = b2%det(1:N_int,2,i2)
i2=i2+1
else if (i2 > b2%cur) then
val(i) = b1%val(i1)
detmp(1:N_int,1,i) = b1%det(1:N_int,1,i1)
detmp(1:N_int,2,i) = b1%det(1:N_int,2,i1)
i1=i1+1
else
if (b1%val(i1) <= b2%val(i2)) then
val(i) = b1%val(i1)
detmp(1:N_int,1,i) = b1%det(1:N_int,1,i1)
detmp(1:N_int,2,i) = b1%det(1:N_int,2,i1)
i1=i1+1
else
val(i) = b2%val(i2)
detmp(1:N_int,1,i) = b2%det(1:N_int,1,i2)
detmp(1:N_int,2,i) = b2%det(1:N_int,2,i2)
i2=i2+1
endif
endif
enddo
deallocate(b2%det, b2%val)
2017-11-26 11:09:03 +01:00
do i=nmwen+1,b2%N
val(i) = 0.d0
detmp(1:N_int,1:2,i) = 0_bit_kind
enddo
2017-05-05 11:32:17 +02:00
b2%det => detmp
b2%val => val
b2%mini = min(b2%mini,b2%val(b2%N))
b2%cur = nmwen
end
2016-09-01 14:43:13 +02:00
subroutine sort_selection_buffer(b)
use selection_types
implicit none
type(selection_buffer), intent(inout) :: b
integer, allocatable :: iorder(:)
2017-03-30 01:13:27 +02:00
integer(bit_kind), pointer :: detmp(:,:,:)
2016-09-01 14:43:13 +02:00
integer :: i, nmwen
logical, external :: detEq
if (b%N == 0 .or. b%cur == 0) return
2016-09-01 14:43:13 +02:00
nmwen = min(b%N, b%cur)
2017-05-05 10:21:31 +02:00
allocate(iorder(b%cur), detmp(N_int, 2, size(b%det,3)))
2016-09-01 14:43:13 +02:00
do i=1,b%cur
iorder(i) = i
end do
2017-05-05 10:21:31 +02:00
call dsort(b%val, 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
end do
2017-05-05 11:32:17 +02:00
deallocate(b%det,iorder)
2017-03-30 01:13:27 +02:00
b%det => detmp
2017-05-05 10:21:31 +02:00
b%mini = min(b%mini,b%val(b%N))
2016-09-01 14:43:13 +02:00
b%cur = nmwen
end subroutine