10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-09-27 03:51:01 +02:00

Fixed tests

This commit is contained in:
Anthony Scemama 2018-12-19 14:33:42 +01:00
parent d1889024d0
commit 1e2255ac88
8 changed files with 37 additions and 23 deletions

View File

@ -12,7 +12,7 @@ Please source the quantum_package.rc file."
;;
let bit_kind_size = lazy (
let filename = root^"/src/Bitmask/bitmasks_module.f90" in
let filename = root^"/src/bitmask/bitmasks_module.f90" in
if not (Sys.file_exists_exn filename) then
raise (Failure ("File "^filename^" not found"));

View File

@ -488,7 +488,7 @@ subroutine save_wavefunction
endif
if (mpi_master) then
call save_wavefunction_general(N_det,N_states,psi_det_sorted,size(psi_coef_sorted,1),psi_coef_sorted)
endif
endif
end
@ -526,7 +526,6 @@ subroutine save_wavefunction_general(ndet,nstates,psidet,dim_psicoef,psicoef)
call ezfio_set_determinants_mo_label(mo_label)
allocate (psi_det_save(N_int,2,ndet))
!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(i,j,k) SHARED(psi_det_save,psidet,ndet,N_int,accu_norm)
do i=1,ndet
do j=1,2
do k=1,N_int
@ -534,7 +533,6 @@ subroutine save_wavefunction_general(ndet,nstates,psidet,dim_psicoef,psicoef)
enddo
enddo
enddo
!$OMP END PARALLEL DO
call ezfio_set_determinants_psi_det(psi_det_save)
deallocate (psi_det_save)

View File

@ -23,9 +23,9 @@ subroutine write_time(iunit)
call print_memory_usage()
call cpu_time(ct)
call wall_time(wt)
write(6,'(A,F15.6,A,F15.6,A)') &
write(6,'(A,F14.6,A,F14.6,A)') &
'.. >>>>> [ WALL TIME: ', wt-output_wall_time_0, &
' s ] [ CPU TIME: ', ct-output_cpu_time_0, ' s ] <<<<< ..'
' s ] [ CPU TIME: ', ct-output_cpu_time_0, ' s ] <<<<< ..'
write(6,*)
end

View File

@ -315,6 +315,7 @@ subroutine select_singles_and_doubles(i_generator,hole_mask,particle_mask,fock_d
integer(bit_kind), allocatable :: minilist(:, :, :), fullminilist(:, :, :)
logical, allocatable :: banned(:,:,:), bannedOrb(:,:)
double precision, allocatable :: mat(:,:,:)
logical :: monoAdo, monoBdo
@ -664,6 +665,9 @@ subroutine fill_buffer_double(i_generator, sp, h1, h2, bannedOrb, banned, fock_d
double precision :: e_pert, delta_E, val, Hii, sum_e_pert, tmp, alpha_h_psi, coef
double precision, external :: diag_H_mat_elem_fock
double precision :: E_shift
! double precision, allocatable :: mat_inv(:,:,:)
!
! allocate(mat_inv(N_states,mo_tot_num,mo_tot_num))
logical, external :: detEq
@ -687,6 +691,7 @@ subroutine fill_buffer_double(i_generator, sp, h1, h2, bannedOrb, banned, fock_d
if(bannedOrb(p1, s1)) cycle
ib = 1
if(sp /= 3) ib = p1+1
! mat_inv(1:N_states,p1,ib:mo_tot_num) = 1.d0/mat(1:N_states,p1,ib:mo_tot_num)
do p2=ib,mo_tot_num
if(bannedOrb(p2, s2)) cycle
if(banned(p1,p2)) cycle
@ -707,7 +712,8 @@ subroutine fill_buffer_double(i_generator, sp, h1, h2, bannedOrb, banned, fock_d
tmp = -tmp
endif
e_pert = 0.5d0 * (tmp - delta_E)
coef = alpha_h_psi / delta_E
coef = e_pert / alpha_h_psi
! coef = e_pert * mat_inv(istate,p1,p2)
pt2(istate) = pt2(istate) + e_pert
variance(istate) = variance(istate) + alpha_h_psi * alpha_h_psi
norm(istate) = norm(istate) + coef * coef
@ -1194,18 +1200,17 @@ subroutine get_d0(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)
logical :: ok
integer :: bant
double precision, allocatable :: hij_cache(:,:)
allocate (hij_cache(mo_tot_num,2))
double precision, allocatable :: hij_cache1(:), hij_cache2(:)
allocate (hij_cache1(mo_tot_num),hij_cache2(mo_tot_num))
bant = 1
if(sp == 3) then ! AB
h1 = p(1,1)
h2 = p(1,2)
do p2=1, mo_tot_num
if(bannedOrb(p2,2)) cycle
call get_mo_bielec_integrals(p2,h1,h2,mo_tot_num,hij_cache(1,1),mo_integrals_map)
call get_mo_bielec_integrals(p2,h1,h2,mo_tot_num,hij_cache1,mo_integrals_map)
do p1=1, mo_tot_num
if(bannedOrb(p1, 1)) cycle
if(banned(p1, p2, bant)) cycle ! rentable?
@ -1214,11 +1219,21 @@ subroutine get_d0(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)
call i_h_j(gen, det, N_int, hij)
else
phase = get_phase_bi(phasemask, 1, 2, h1, p1, h2, p2, N_int)
hij = hij_cache(p1,1) * phase
hij = hij_cache1(p1) * phase
end if
do k=1,N_states
mat(k, p1, p2) = mat(k, p1, p2) + coefs(k) * hij ! HOTSPOT
enddo
! if( (bannedOrb(p1, 1)).or.(banned(p1, p2, bant)) ) then
! hij = 0.d0
! else if(p1 /= h1 .and. p2 /= h2) then
! hij = hij_cache1(p1) * get_phase_bi(phasemask, 1, 2, h1, p1, h2, p2, N_int)
! else
! call apply_particles(mask, 1,p1,2,p2, det, ok, N_int)
! call i_h_j(gen, det, N_int, hij)
! end if
if (hij /= 0.d0) then
do k=1,N_states
mat(k, p1, p2) = mat(k, p1, p2) + coefs(k) * hij ! HOTSPOT
enddo
endif
end do
end do
@ -1227,8 +1242,8 @@ subroutine get_d0(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)
p2 = p(2,sp)
do puti=1, mo_tot_num
if(bannedOrb(puti, sp)) cycle
call get_mo_bielec_integrals(puti,p2,p1,mo_tot_num,hij_cache(1,1),mo_integrals_map)
call get_mo_bielec_integrals(puti,p1,p2,mo_tot_num,hij_cache(1,2),mo_integrals_map)
call get_mo_bielec_integrals(puti,p2,p1,mo_tot_num,hij_cache1,mo_integrals_map)
call get_mo_bielec_integrals(puti,p1,p2,mo_tot_num,hij_cache2,mo_integrals_map)
do putj=puti+1, mo_tot_num
if(bannedOrb(putj, sp)) cycle
if(banned(puti, putj, bant)) cycle ! rentable?
@ -1241,7 +1256,7 @@ subroutine get_d0(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)
enddo
endif
else
hij = hij_cache(putj,1) - hij_cache(putj,2)
hij = hij_cache1(putj) - hij_cache2(putj)
if (hij /= 0.d0) then
hij = hij * get_phase_bi(phasemask, sp, sp, puti, p1 , putj, p2, N_int)
do k=1,N_states
@ -1252,7 +1267,8 @@ subroutine get_d0(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)
end do
end do
end if
deallocate(hij_cache)
deallocate(hij_cache1,hij_cache2)
end

Binary file not shown.

View File

@ -4,7 +4,7 @@ source $QP_ROOT/tests/bats/common.bats.sh
function run_FCI() {
thresh=5.e-5
test_exe FCI || skip
test_exe fci || skip
qp_edit -c $1
ezfio set_file $1
ezfio set perturbation do_pt2 True

View File

@ -11,7 +11,7 @@ function run_init() {
function run_HF() {
thresh=1.e-8
test_exe SCF || skip
test_exe scf || skip
qp_edit -c $1
ezfio set_file $1
ezfio set hartree_fock thresh_scf 1.e-10

View File

@ -11,7 +11,7 @@ function run_init() {
function run_HF() {
thresh=1.e-7
test_exe SCF || skip
test_exe scf || skip
qp_edit -c $1
ezfio set_file $1
ezfio set hartree_fock thresh_scf 2.e-8
@ -23,7 +23,7 @@ function run_HF() {
function run_FCI() {
thresh=5.e-5
test_exe FCI || skip
test_exe fci || skip
qp_edit -c $1
ezfio set_file $1
ezfio set perturbation do_pt2 True