10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-26 07:02:14 +02:00

Accelerated selection

This commit is contained in:
Anthony Scemama 2017-04-19 15:31:12 +02:00
parent 69a747fde0
commit 1ac36ab762
3 changed files with 50 additions and 30 deletions

View File

@ -45,7 +45,7 @@ subroutine run_selection_slave(thread,iproc,energy)
if(buf%N == 0) then
! Only first time
call create_selection_buffer(N, N*2, buf)
call create_selection_buffer(N, N*3, buf2)
call create_selection_buffer(N, N*2, buf2)
else
if(N /= buf%N) stop "N changed... wtf man??"
end if
@ -62,7 +62,6 @@ subroutine run_selection_slave(thread,iproc,energy)
do i=1,buf%cur
call add_to_selection_buffer(buf2, buf%det(1,1,i), buf%val(i))
enddo
call sort_selection_buffer(buf2)
buf%mini = buf2%mini
pt2 = 0d0
buf%cur = 0

View File

@ -41,7 +41,6 @@ subroutine sort_selection_buffer(b)
implicit none
type(selection_buffer), intent(inout) :: b
double precision, allocatable:: absval(:)
integer, allocatable :: iorder(:)
double precision, pointer :: vals(:)
integer(bit_kind), pointer :: detmp(:,:,:)
@ -49,29 +48,23 @@ subroutine sort_selection_buffer(b)
logical, external :: detEq
nmwen = min(b%N, b%cur)
allocate(iorder(b%cur), detmp(N_int, 2, size(b%det,3)), absval(b%cur), vals(size(b%val)))
absval = -dabs(b%val(:b%cur))
allocate(iorder(b%cur), detmp(N_int, 2, size(b%det,3)), vals(size(b%val)))
do i=1,b%cur
iorder(i) = i
end do
! Optimal for almost sorted data
! call sorted_dnumber(absval, b%cur, i)
! if (b%cur/i >
! call insertion_dsort(absval, iorder, b%cur)
call dsort(absval, iorder, b%cur)
call dsort(b%val, iorder, b%cur)
do i=1, nmwen
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))
vals(i) = b%val(iorder(i))
end do
do i=nmwen+1, size(vals)
vals(i) = 0.d0
enddo
if (nmwen < b%N) then
vals(nmwen+1) = 0.d0
endif
deallocate(b%det, b%val)
b%det => detmp
b%val => vals
b%mini = max(b%mini,dabs(b%val(b%N)))
b%mini = min(b%mini,b%val(1))
b%cur = nmwen
end subroutine

View File

@ -12,23 +12,20 @@ BEGIN_TEMPLATE
$type :: xtmp
integer :: i, i0, j, jmax
do i=1,isize
do i=2,isize
xtmp = x(i)
i0 = iorder(i)
j = i-1
do j=i-1,1,-1
if ( x(j) > xtmp ) then
x(j+1) = x(j)
iorder(j+1) = iorder(j)
else
exit
endif
do j = i-1,1,-1
if ((x(j) <= xtmp)) exit
x(j+1) = x(j)
iorder(j+1) = iorder(j)
enddo
x(j+1) = xtmp
iorder(j+1) = i0
enddo
end subroutine insertion_$Xsort
subroutine heap_$Xsort(x,iorder,isize)
implicit none
BEGIN_DOC
@ -179,7 +176,7 @@ BEGIN_TEMPLATE
endif
do i=2,isize
if (x(i-1) >= x(i)) then
if (x(i-1) <= x(i)) then
n=n+1
endif
enddo
@ -194,6 +191,31 @@ SUBST [ X, type ]
i2 ; integer*2 ;;
END_TEMPLATE
!BEGIN_TEMPLATE
! subroutine $Xsort(x,iorder,isize)
! implicit none
! BEGIN_DOC
! ! Sort array x(isize).
! ! iorder in input should be (1,2,3,...,isize), and in output
! ! contains the new order of the elements.
! END_DOC
! integer,intent(in) :: isize
! $type,intent(inout) :: x(isize)
! integer,intent(inout) :: iorder(isize)
! integer :: n
! call sorted_$Xnumber(x,isize,n)
! if ( isize-n < 1000) then
! call insertion_$Xsort(x,iorder,isize)
! else
! call heap_$Xsort(x,iorder,isize)
! endif
! end subroutine $Xsort
!
!SUBST [ X, type ]
! ; real ;;
! d ; double precision ;;
!END_TEMPLATE
BEGIN_TEMPLATE
subroutine $Xsort(x,iorder,isize)
implicit none
@ -207,16 +229,19 @@ BEGIN_TEMPLATE
integer,intent(inout) :: iorder(isize)
integer :: n
call sorted_$Xnumber(x,isize,n)
if ( isize-n < 1000) then
if (isize == n) then
return
endif
if ( isize < 512+n) then
call insertion_$Xsort(x,iorder,isize)
else
call heap_$Xsort(x,iorder,isize)
call $Yradix_sort(x,iorder,isize,-1)
endif
end subroutine $Xsort
SUBST [ X, type ]
; real ;;
d ; double precision ;;
SUBST [ X, type, Y ]
; real ; i ;;
d ; double precision ; i8 ;;
END_TEMPLATE
BEGIN_TEMPLATE
@ -422,6 +447,9 @@ BEGIN_TEMPLATE
i0 = 0_$int_type
i4 = maxval(x)
if (i4 == 0_$type) then
return
endif
iradix_new = $integer_size-1-leadz(i4)
mask = ibset(0_$type,iradix_new)