10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-22 20:35:19 +01:00

Fixed make_s2_eigenfunction

This commit is contained in:
Anthony Scemama 2016-11-08 10:24:36 +01:00
parent af2780860e
commit 7ac373c1b3
3 changed files with 44 additions and 23 deletions

View File

@ -122,13 +122,15 @@ subroutine ZMQ_selection(N_in, pt2)
double precision, intent(out) :: pt2(N_states)
N = max(N_in,1)
provide nproc
provide ci_electronic_energy
call new_parallel_job(zmq_to_qp_run_socket,"selection")
call zmq_put_psi(zmq_to_qp_run_socket,1,ci_electronic_energy,size(ci_electronic_energy))
call zmq_set_running(zmq_to_qp_run_socket)
call create_selection_buffer(N, N*2, b)
if (.True.) then
N = max(N_in,1)
provide nproc
provide ci_electronic_energy
call new_parallel_job(zmq_to_qp_run_socket,"selection")
call zmq_put_psi(zmq_to_qp_run_socket,1,ci_electronic_energy,size(ci_electronic_energy))
call zmq_set_running(zmq_to_qp_run_socket)
call create_selection_buffer(N, N*2, b)
endif
integer :: i_generator, i_generator_start, i_generator_max, step
! step = int(max(1.,10*elec_num/mo_tot_num)
@ -154,6 +156,7 @@ subroutine ZMQ_selection(N_in, pt2)
if (N_in > 0) then
call fill_H_apply_buffer_no_selection(b%cur,b%det,N_int,0) !!! PAS DE ROBIN
call copy_H_apply_buffer_to_wf()
call make_s2_eigenfunction
endif
end subroutine

View File

@ -259,6 +259,7 @@ subroutine remove_duplicates_in_psi_det(found_duplicates)
N_det = k
call write_bool(output_determinants,found_duplicates,'Found duplicate determinants')
SOFT_TOUCH N_det psi_det psi_coef
stop 'duplicates in psi_det'
endif
deallocate (duplicate,bit_tmp)
end

View File

@ -35,7 +35,8 @@ subroutine occ_pattern_to_dets_size(o,sze,n_alpha,Nint)
bmax += popcnt( o(k,1) )
amax -= popcnt( o(k,2) )
enddo
sze = 2*int( min(binom_func(bmax, amax), 1.d8) )
sze = int( min(binom_func(bmax, amax), 1.d8) )
sze = sze*sze
end
@ -123,8 +124,8 @@ end
implicit none
BEGIN_DOC
! array of the occ_pattern present in the wf
! psi_occ_pattern(:,1,j) = jth occ_pattern of the wave function : represent all the single occupation
! psi_occ_pattern(:,2,j) = jth occ_pattern of the wave function : represent all the double occupation
! psi_occ_pattern(:,1,j) = jth occ_pattern of the wave function : represent all the single occupations
! psi_occ_pattern(:,2,j) = jth occ_pattern of the wave function : represent all the double occupations
END_DOC
integer :: i,j,k
@ -144,7 +145,7 @@ end
logical,allocatable :: duplicate(:)
allocate ( iorder(N_det), duplicate(N_det), bit_tmp(N_det), tmp_array(N_int,2,psi_det_size) )
allocate ( iorder(N_det), duplicate(N_det), bit_tmp(N_det), tmp_array(N_int,2,N_det) )
do i=1,N_det
iorder(i) = i
@ -161,18 +162,16 @@ end
duplicate(i) = .False.
enddo
i=1
integer (bit_kind) :: occ_pattern_tmp
do i=1,N_det
duplicate(i) = .False.
enddo
! Find duplicates
do i=1,N_det-1
if (duplicate(i)) then
cycle
endif
j = i+1
do while (bit_tmp(j)==bit_tmp(i))
if (j>N_det) then
exit
endif
if (duplicate(j)) then
j+=1
cycle
@ -186,12 +185,10 @@ end
endif
enddo
j+=1
if (j>=N_det) then
exit
endif
enddo
enddo
! Copy filtered result
N_occ_pattern=0
do i=1,N_det
if (duplicate(i)) then
@ -204,6 +201,28 @@ end
enddo
enddo
!- Check
! do i=1,N_occ_pattern
! do j=i+1,N_occ_pattern
! duplicate(1) = .True.
! do k=1,N_int
! if (psi_occ_pattern(k,1,i) /= psi_occ_pattern(k,1,j)) then
! duplicate(1) = .False.
! exit
! endif
! if (psi_occ_pattern(k,2,i) /= psi_occ_pattern(k,2,j)) then
! duplicate(1) = .False.
! exit
! endif
! enddo
! if (duplicate(1)) then
! call debug_det(psi_occ_pattern(1,1,i),N_int)
! call debug_det(psi_occ_pattern(1,1,j),N_int)
! stop 'DUPLICATE'
! endif
! enddo
! enddo
!-
deallocate(iorder,duplicate,bit_tmp,tmp_array)
END_PROVIDER
@ -217,9 +236,6 @@ subroutine make_s2_eigenfunction
integer, parameter :: bufsze = 1000
logical, external :: is_in_wavefunction
return
stop 'make_s2_eigenfunction has a bug. It should not be used!!!'
allocate (d(N_int,2,1), det_buffer(N_int,2,bufsze) )
smax = 1
N_det_new = 0
@ -258,6 +274,7 @@ subroutine make_s2_eigenfunction
call copy_H_apply_buffer_to_wf
SOFT_TOUCH N_det psi_coef psi_det
print *, 'Added determinants for S^2'
! call remove_duplicates_in_psi_det
end