mirror of
https://github.com/LCPQ/quantum_package
synced 2024-12-22 12:23:48 +01:00
Fixed CASSD
This commit is contained in:
parent
11aeaa91c7
commit
e0183e998c
@ -41,8 +41,8 @@ subroutine run_selection_slave(thread,iproc,energy)
|
||||
if (done) then
|
||||
ctask = ctask - 1
|
||||
else
|
||||
integer :: i_generator, i_generator_start, i_generator_max, step, N
|
||||
read (task,*) i_generator_start, i_generator_max, step, N
|
||||
integer :: i_generator, N
|
||||
read (task,*) i_generator, N
|
||||
if(buf%N == 0) then
|
||||
! Only first time
|
||||
call create_selection_buffer(N, N*2, buf)
|
||||
@ -50,9 +50,7 @@ subroutine run_selection_slave(thread,iproc,energy)
|
||||
else
|
||||
if(N /= buf%N) stop "N changed... wtf man??"
|
||||
end if
|
||||
do i_generator=i_generator_start,i_generator_max,step
|
||||
call select_connected(i_generator,energy,pt2,buf)
|
||||
enddo
|
||||
call select_connected(i_generator,energy,pt2,buf)
|
||||
endif
|
||||
|
||||
if(done .or. ctask == size(task_id)) then
|
||||
|
@ -1215,34 +1215,41 @@ subroutine ZMQ_selection(N_in, pt2)
|
||||
|
||||
implicit none
|
||||
|
||||
character*(512) :: task
|
||||
integer(ZMQ_PTR) :: zmq_to_qp_run_socket
|
||||
integer, intent(in) :: N_in
|
||||
type(selection_buffer) :: b
|
||||
integer :: i, N
|
||||
integer, external :: omp_get_thread_num
|
||||
double precision, intent(out) :: pt2(N_states)
|
||||
integer, parameter :: maxtasks=10000
|
||||
|
||||
|
||||
N = max(N_in,1)
|
||||
if (.True.) then
|
||||
PROVIDE pt2_e0_denominator
|
||||
N = max(N_in,1)
|
||||
provide nproc
|
||||
call new_parallel_job(zmq_to_qp_run_socket,"selection")
|
||||
call zmq_put_psi(zmq_to_qp_run_socket,1,pt2_e0_denominator,size(pt2_e0_denominator))
|
||||
call create_selection_buffer(N, N*2, b)
|
||||
endif
|
||||
|
||||
integer :: i_generator, i_generator_start, i_generator_max, step
|
||||
character*(20*maxtasks) :: task
|
||||
task = ' '
|
||||
|
||||
step = int(5000000.d0 / dble(N_int * N_states * elec_num * elec_num * mo_tot_num * mo_tot_num ))
|
||||
step = max(1,step)
|
||||
do i= 1, N_det_generators,step
|
||||
i_generator_start = i
|
||||
i_generator_max = min(i+step-1,N_det_generators)
|
||||
write(task,*) i_generator_start, i_generator_max, 1, N
|
||||
integer :: k
|
||||
k=0
|
||||
do i= 1, N_det_generators
|
||||
k = k+1
|
||||
write(task(20*(k-1)+1:20*k),'(I9,1X,I9,''|'')') i, N
|
||||
k = k+20
|
||||
if (k>20*maxtasks) then
|
||||
k=0
|
||||
call add_task_to_taskserver(zmq_to_qp_run_socket,task)
|
||||
endif
|
||||
enddo
|
||||
if (k > 0) then
|
||||
call add_task_to_taskserver(zmq_to_qp_run_socket,task)
|
||||
end do
|
||||
endif
|
||||
call zmq_set_running(zmq_to_qp_run_socket)
|
||||
|
||||
!$OMP PARALLEL DEFAULT(shared) SHARED(b, pt2) PRIVATE(i) NUM_THREADS(nproc+1)
|
||||
@ -1250,7 +1257,6 @@ subroutine ZMQ_selection(N_in, pt2)
|
||||
if (i==0) then
|
||||
call selection_collector(b, pt2)
|
||||
else
|
||||
call sleep(1)
|
||||
call selection_slave_inproc(i)
|
||||
endif
|
||||
!$OMP END PARALLEL
|
||||
@ -1261,6 +1267,7 @@ subroutine ZMQ_selection(N_in, pt2)
|
||||
if (s2_eig) then
|
||||
call make_s2_eigenfunction
|
||||
endif
|
||||
call save_wavefunction
|
||||
endif
|
||||
end subroutine
|
||||
|
||||
|
@ -24,7 +24,6 @@ subroutine ZMQ_selection(N_in, pt2)
|
||||
call create_selection_buffer(N, N*2, b)
|
||||
endif
|
||||
|
||||
! Ugly, but variable-length strings don't work as expected with gfortran < 4.8 :-(
|
||||
character*(20*maxtasks) :: task
|
||||
task = ' '
|
||||
|
||||
|
@ -44,8 +44,8 @@ subroutine bielec_integrals_index_reverse(i,j,k,l,i1)
|
||||
l(1) = ceiling(0.5d0*(dsqrt(8.d0*dble(i2)+1.d0)-1.d0))
|
||||
i3 = i1 - ishft(i2*i2-i2,-1)
|
||||
k(1) = ceiling(0.5d0*(dsqrt(8.d0*dble(i3)+1.d0)-1.d0))
|
||||
j(1) = i2 - ishft(l(1)*l(1)-l(1),-1)
|
||||
i(1) = i3 - ishft(k(1)*k(1)-k(1),-1)
|
||||
j(1) = int(i2 - ishft(l(1)*l(1)-l(1),-1),4)
|
||||
i(1) = int(i3 - ishft(k(1)*k(1)-k(1),-1),4)
|
||||
|
||||
!ijkl
|
||||
i(2) = i(1) !ilkj
|
||||
|
@ -26,7 +26,7 @@ subroutine svd(A,LDA,U,LDU,D,Vt,LDVt,m,n)
|
||||
lwork = -1
|
||||
call dgesvd('A','A', m, n, A_tmp, LDA, &
|
||||
D, U, LDU, Vt, LDVt, work, lwork, info)
|
||||
lwork = work(1)
|
||||
lwork = int(work(1))
|
||||
deallocate(work)
|
||||
|
||||
allocate(work(lwork))
|
||||
@ -149,7 +149,7 @@ subroutine ortho_qr(A,LDA,m,n)
|
||||
allocate (jpvt(n), tau(n), work(1))
|
||||
LWORK=-1
|
||||
call dgeqrf( m, n, A, LDA, TAU, WORK, LWORK, INFO )
|
||||
LWORK=2*WORK(1)
|
||||
LWORK=2*int(WORK(1))
|
||||
deallocate(WORK)
|
||||
allocate(WORK(LWORK))
|
||||
call dgeqrf(m, n, A, LDA, TAU, WORK, LWORK, INFO )
|
||||
@ -293,7 +293,7 @@ subroutine get_pseudo_inverse(A,m,n,C,LDA)
|
||||
print *, info, ': SVD failed'
|
||||
stop
|
||||
endif
|
||||
lwork = work(1)
|
||||
lwork = int(work(1))
|
||||
deallocate(work)
|
||||
allocate(work(lwork))
|
||||
call dgesvd('S','A', m, n, A_tmp, m,D,U,m,Vt,n,work,lwork,info)
|
||||
|
@ -105,7 +105,7 @@ subroutine map_load_from_disk(filename,map)
|
||||
map % map(i) % value => map % consolidated_value ( map % consolidated_idx (i+1) :)
|
||||
map % map(i) % key => map % consolidated_key ( map % consolidated_idx (i+1) :)
|
||||
map % map(i) % sorted = .True.
|
||||
n_elements = map % consolidated_idx (i+2) - k
|
||||
n_elements = int( map % consolidated_idx (i+2) - k, 4)
|
||||
k = map % consolidated_idx (i+2)
|
||||
map % map(i) % map_size = n_elements
|
||||
map % map(i) % n_elements = n_elements
|
||||
|
@ -53,17 +53,17 @@ module map_module
|
||||
end module map_module
|
||||
|
||||
|
||||
real function map_mb(map)
|
||||
double precision function map_mb(map)
|
||||
use map_module
|
||||
use omp_lib
|
||||
implicit none
|
||||
type (map_type), intent(in) :: map
|
||||
integer(map_size_kind) :: i
|
||||
|
||||
map_mb = 8+map_size_kind+map_size_kind+omp_lock_kind+4
|
||||
map_mb = dble(8+map_size_kind+map_size_kind+omp_lock_kind+4)
|
||||
do i=0,map%map_size
|
||||
map_mb = map_mb + map%map(i)%map_size*(cache_key_kind+integral_kind) +&
|
||||
8+8+4+cache_map_size_kind+cache_map_size_kind+omp_lock_kind
|
||||
map_mb = map_mb + dble(map%map(i)%map_size*(cache_key_kind+integral_kind) +&
|
||||
8+8+4+cache_map_size_kind+cache_map_size_kind+omp_lock_kind)
|
||||
enddo
|
||||
map_mb = map_mb / (1024.d0*1024.d0)
|
||||
end
|
||||
@ -406,7 +406,7 @@ subroutine map_update(map, key, value, sze, thr)
|
||||
call cache_map_reallocate(local_map, local_map%n_elements + local_map%n_elements)
|
||||
call cache_map_shrink(local_map,thr)
|
||||
endif
|
||||
cache_key = iand(key(i),map_mask)
|
||||
cache_key = int(iand(key(i),map_mask),2)
|
||||
local_map%n_elements = local_map%n_elements + 1
|
||||
local_map%value(local_map%n_elements) = value(i)
|
||||
local_map%key(local_map%n_elements) = cache_key
|
||||
@ -464,7 +464,7 @@ subroutine map_append(map, key, value, sze)
|
||||
if (n_elements == map%map(idx_cache)%map_size) then
|
||||
call cache_map_reallocate(map%map(idx_cache), n_elements+ ishft(n_elements,-1))
|
||||
endif
|
||||
cache_key = iand(key(i),map_mask)
|
||||
cache_key = int(iand(key(i),map_mask),2)
|
||||
map%map(idx_cache)%value(n_elements) = value(i)
|
||||
map%map(idx_cache)%key(n_elements) = cache_key
|
||||
map%map(idx_cache)%n_elements = n_elements
|
||||
@ -615,7 +615,7 @@ subroutine search_key_big_interval(key,X,sze,idx,ibegin_in,iend_in)
|
||||
idx = -1
|
||||
return
|
||||
endif
|
||||
cache_key = iand(key,map_mask)
|
||||
cache_key = int(iand(key,map_mask),2)
|
||||
ibegin = min(ibegin_in,sze)
|
||||
iend = min(iend_in,sze)
|
||||
if ((cache_key > X(ibegin)) .and. (cache_key < X(iend))) then
|
||||
@ -723,7 +723,7 @@ subroutine search_key_value_big_interval(key,value,X,Y,sze,idx,ibegin_in,iend_in
|
||||
value = 0.d0
|
||||
return
|
||||
endif
|
||||
cache_key = iand(key,map_mask)
|
||||
cache_key = int(iand(key,map_mask),2)
|
||||
ibegin = min(ibegin_in,sze)
|
||||
iend = min(iend_in,sze)
|
||||
if ((cache_key > X(ibegin)) .and. (cache_key < X(iend))) then
|
||||
|
@ -292,18 +292,17 @@ BEGIN_TEMPLATE
|
||||
! contains the new order of the elements.
|
||||
! iradix should be -1 in input.
|
||||
END_DOC
|
||||
$int_type, intent(in) :: isize
|
||||
$int_type, intent(inout) :: iorder(isize)
|
||||
$type, intent(inout) :: x(isize)
|
||||
integer*$int_type, intent(in) :: isize
|
||||
integer*$int_type, intent(inout) :: iorder(isize)
|
||||
integer*$type, intent(inout) :: x(isize)
|
||||
integer, intent(in) :: iradix
|
||||
integer :: iradix_new
|
||||
$type, allocatable :: x2(:), x1(:)
|
||||
$type :: i4
|
||||
$int_type, allocatable :: iorder1(:),iorder2(:)
|
||||
$int_type :: i0, i1, i2, i3, i
|
||||
integer*$type, allocatable :: x2(:), x1(:)
|
||||
integer*$type :: i4
|
||||
integer*$int_type, allocatable :: iorder1(:),iorder2(:)
|
||||
integer*$int_type :: i0, i1, i2, i3, i
|
||||
integer, parameter :: integer_size=$octets
|
||||
$type, parameter :: zero=$zero
|
||||
$type :: mask
|
||||
integer*$type :: mask
|
||||
integer :: nthreads, omp_get_num_threads
|
||||
!DIR$ ATTRIBUTES ALIGN : 128 :: iorder1,iorder2, x2, x1
|
||||
|
||||
@ -311,16 +310,16 @@ BEGIN_TEMPLATE
|
||||
|
||||
! Find most significant bit
|
||||
|
||||
i0 = 0_8
|
||||
i4 = -1_8
|
||||
i0 = 0_$int_type
|
||||
i4 = -1_$type
|
||||
|
||||
do i=1,isize
|
||||
i4 = max(i4,x(i))
|
||||
enddo
|
||||
i3 = i4 ! Type conversion
|
||||
i3 = int(i4,$int_type)
|
||||
|
||||
iradix_new = integer_size-1-leadz(i3)
|
||||
mask = ibset(zero,iradix_new)
|
||||
mask = ibset(0_$type,iradix_new)
|
||||
nthreads = 1
|
||||
! nthreads = 1+ishft(omp_get_num_threads(),-1)
|
||||
|
||||
@ -331,22 +330,22 @@ BEGIN_TEMPLATE
|
||||
stop
|
||||
endif
|
||||
|
||||
i1=1_8
|
||||
i2=1_8
|
||||
i1=1_$int_type
|
||||
i2=1_$int_type
|
||||
|
||||
do i=1,isize
|
||||
if (iand(mask,x(i)) == zero) then
|
||||
if (iand(mask,x(i)) == 0_$type) then
|
||||
iorder1(i1) = iorder(i)
|
||||
x1(i1) = x(i)
|
||||
i1 = i1+1_8
|
||||
i1 = i1+1_$int_type
|
||||
else
|
||||
iorder2(i2) = iorder(i)
|
||||
x2(i2) = x(i)
|
||||
i2 = i2+1_8
|
||||
i2 = i2+1_$int_type
|
||||
endif
|
||||
enddo
|
||||
i1=i1-1_8
|
||||
i2=i2-1_8
|
||||
i1=i1-1_$int_type
|
||||
i2=i2-1_$int_type
|
||||
|
||||
do i=1,i1
|
||||
iorder(i0+i) = iorder1(i)
|
||||
@ -399,12 +398,12 @@ BEGIN_TEMPLATE
|
||||
endif
|
||||
|
||||
|
||||
mask = ibset(zero,iradix)
|
||||
mask = ibset(0_$type,iradix)
|
||||
i0=1
|
||||
i1=1
|
||||
|
||||
do i=1,isize
|
||||
if (iand(mask,x(i)) == zero) then
|
||||
if (iand(mask,x(i)) == 0_$type) then
|
||||
iorder(i0) = iorder(i)
|
||||
x(i0) = x(i)
|
||||
i0 = i0+1
|
||||
@ -443,12 +442,12 @@ BEGIN_TEMPLATE
|
||||
|
||||
end
|
||||
|
||||
SUBST [ X, type, octets, is_big, big, int_type, zero ]
|
||||
i ; integer ; 32 ; .False. ; ; integer ; 0;;
|
||||
i8 ; integer*8 ; 32 ; .False. ; ; integer ; 0_8;;
|
||||
i2 ; integer*2 ; 32 ; .False. ; ; integer ; 0;;
|
||||
i ; integer ; 64 ; .True. ; _big ; integer*8 ; 0 ;;
|
||||
i8 ; integer*8 ; 64 ; .True. ; _big ; integer*8 ; 0_8 ;;
|
||||
SUBST [ X, type, octets, is_big, big, int_type ]
|
||||
i ; 4 ; 32 ; .False. ; ; 4 ;;
|
||||
i8 ; 8 ; 32 ; .False. ; ; 4 ;;
|
||||
i2 ; 2 ; 32 ; .False. ; ; 4 ;;
|
||||
i ; 4 ; 64 ; .True. ; _big ; 8 ;;
|
||||
i8 ; 8 ; 64 ; .True. ; _big ; 8 ;;
|
||||
END_TEMPLATE
|
||||
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user