mirror of
https://github.com/LCPQ/quantum_package
synced 2024-11-03 20:54:00 +01:00
commit
edcfa8be85
6
plugins/QMC/EZFIO.cfg
Normal file
6
plugins/QMC/EZFIO.cfg
Normal file
@ -0,0 +1,6 @@
|
|||||||
|
[ci_threshold]
|
||||||
|
type: Threshold
|
||||||
|
doc: Threshold on the CI coefficients as in QMCChem
|
||||||
|
interface: ezfio,provider,ocaml
|
||||||
|
default: 1.e-8
|
||||||
|
|
@ -1 +1 @@
|
|||||||
Determinants Davidson
|
Determinants Davidson Full_CI_ZMQ
|
||||||
|
@ -24,13 +24,13 @@ program save_for_qmc
|
|||||||
)
|
)
|
||||||
iunit = 13
|
iunit = 13
|
||||||
open(unit=iunit,file=trim(ezfio_filename)//'/simulation/e_ref',action='write')
|
open(unit=iunit,file=trim(ezfio_filename)//'/simulation/e_ref',action='write')
|
||||||
call ezfio_has_full_ci_energy_pt2(exists)
|
call ezfio_has_full_ci_zmq_energy_pt2(exists)
|
||||||
if (exists) then
|
if (exists) then
|
||||||
call ezfio_get_full_ci_energy_pt2(e_ref)
|
call ezfio_get_full_ci_zmq_energy_pt2(e_ref)
|
||||||
else
|
else
|
||||||
call ezfio_has_full_ci_energy(exists)
|
call ezfio_has_full_ci_zmq_energy(exists)
|
||||||
if (exists) then
|
if (exists) then
|
||||||
call ezfio_get_full_ci_energy(e_ref)
|
call ezfio_get_full_ci_zmq_energy(e_ref)
|
||||||
else
|
else
|
||||||
call ezfio_has_hartree_fock_energy(exists)
|
call ezfio_has_hartree_fock_energy(exists)
|
||||||
if (exists) then
|
if (exists) then
|
||||||
|
@ -21,6 +21,6 @@ program qmcpack
|
|||||||
enddo
|
enddo
|
||||||
call save_mos
|
call save_mos
|
||||||
call system('rm '//trim(ezfio_filename)//'/mo_basis/ao_md5')
|
call system('rm '//trim(ezfio_filename)//'/mo_basis/ao_md5')
|
||||||
call system('$QP_ROOT/src/qmcpack/qp_convert_qmcpack_to_ezfio.py '//trim(ezfio_filename))
|
call system('$QP_ROOT/src/QMC/qp_convert_qmcpack_to_ezfio.py '//trim(ezfio_filename))
|
||||||
|
|
||||||
end
|
end
|
||||||
|
@ -1,4 +1,10 @@
|
|||||||
program e_curve
|
program truncate
|
||||||
|
read_wf = .True.
|
||||||
|
SOFT_TOUCH read_wf
|
||||||
|
call run
|
||||||
|
end
|
||||||
|
|
||||||
|
subroutine run
|
||||||
use bitmasks
|
use bitmasks
|
||||||
implicit none
|
implicit none
|
||||||
integer :: i,j,k, kk, nab, m, l
|
integer :: i,j,k, kk, nab, m, l
|
||||||
@ -6,24 +12,17 @@ program e_curve
|
|||||||
integer, allocatable :: iorder(:)
|
integer, allocatable :: iorder(:)
|
||||||
double precision , allocatable :: norm_sort(:)
|
double precision , allocatable :: norm_sort(:)
|
||||||
double precision :: e_0(N_states)
|
double precision :: e_0(N_states)
|
||||||
if (.not.read_wf) then
|
|
||||||
stop 'Please set read_wf to true'
|
|
||||||
endif
|
|
||||||
|
|
||||||
PROVIDE mo_bielec_integrals_in_map H_apply_buffer_allocated
|
PROVIDE mo_bielec_integrals_in_map H_apply_buffer_allocated
|
||||||
|
|
||||||
nab = n_det_alpha_unique+n_det_beta_unique
|
nab = n_det_alpha_unique+n_det_beta_unique
|
||||||
allocate ( norm_sort(0:nab), iorder(0:nab) )
|
allocate ( norm_sort(0:nab), iorder(0:nab) )
|
||||||
|
|
||||||
double precision :: thresh
|
|
||||||
integer(bit_kind), allocatable :: det_i(:,:), det_j(:,:)
|
integer(bit_kind), allocatable :: det_i(:,:), det_j(:,:)
|
||||||
double precision, allocatable :: u_t(:,:), v_t(:,:), s_t(:,:)
|
double precision, allocatable :: u_t(:,:), v_t(:,:), s_t(:,:)
|
||||||
double precision, allocatable :: u_0(:,:), v_0(:,:)
|
double precision, allocatable :: u_0(:,:), v_0(:,:)
|
||||||
allocate(u_t(N_states,N_det),v_t(N_states,N_det),s_t(N_states,N_det))
|
allocate(u_t(N_states,N_det),v_t(N_states,N_det),s_t(N_states,N_det))
|
||||||
allocate(u_0(N_states,N_det),v_0(N_states,N_det))
|
allocate(u_0(N_det,N_states),v_0(N_det,N_states))
|
||||||
|
|
||||||
print *, 'Threshold?'
|
|
||||||
read(*,*) thresh
|
|
||||||
|
|
||||||
norm_sort(0) = 0.d0
|
norm_sort(0) = 0.d0
|
||||||
iorder(0) = 0
|
iorder(0) = 0
|
||||||
@ -45,19 +44,23 @@ program e_curve
|
|||||||
do j=0,nab
|
do j=0,nab
|
||||||
i = iorder(j)
|
i = iorder(j)
|
||||||
if (i<0) then
|
if (i<0) then
|
||||||
|
!$OMP PARALLEL DO PRIVATE(k)
|
||||||
do k=1,n_det
|
do k=1,n_det
|
||||||
if (psi_bilinear_matrix_columns(k) == -i) then
|
if (psi_bilinear_matrix_columns(k) == -i) then
|
||||||
psi_bilinear_matrix_values(k,1) = 0.d0
|
psi_bilinear_matrix_values(k,1) = 0.d0
|
||||||
endif
|
endif
|
||||||
enddo
|
enddo
|
||||||
|
!$OMP END PARALLEL DO
|
||||||
else
|
else
|
||||||
|
!$OMP PARALLEL DO PRIVATE(k)
|
||||||
do k=1,n_det
|
do k=1,n_det
|
||||||
if (psi_bilinear_matrix_rows(k) == i) then
|
if (psi_bilinear_matrix_rows(k) == i) then
|
||||||
psi_bilinear_matrix_values(k,1) = 0.d0
|
psi_bilinear_matrix_values(k,1) = 0.d0
|
||||||
endif
|
endif
|
||||||
enddo
|
enddo
|
||||||
|
!$OMP END PARALLEL DO
|
||||||
endif
|
endif
|
||||||
if (thresh > norm_sort(j)) then
|
if (ci_threshold > norm_sort(j)) then
|
||||||
cycle
|
cycle
|
||||||
endif
|
endif
|
||||||
|
|
||||||
@ -70,7 +73,9 @@ program e_curve
|
|||||||
u_t, &
|
u_t, &
|
||||||
size(u_t, 1), &
|
size(u_t, 1), &
|
||||||
N_det, N_states)
|
N_det, N_states)
|
||||||
|
print *, 'Computing H|Psi> ...'
|
||||||
call H_S2_u_0_nstates_openmp_work(v_t,s_t,u_t,N_states,N_det,1,N_det,0,1)
|
call H_S2_u_0_nstates_openmp_work(v_t,s_t,u_t,N_states,N_det,1,N_det,0,1)
|
||||||
|
print *, 'Done'
|
||||||
call dtranspose( &
|
call dtranspose( &
|
||||||
v_t, &
|
v_t, &
|
||||||
size(v_t, 1), &
|
size(v_t, 1), &
|
||||||
@ -96,7 +101,7 @@ program e_curve
|
|||||||
print *, 'Energy', E
|
print *, 'Energy', E
|
||||||
exit
|
exit
|
||||||
enddo
|
enddo
|
||||||
call wf_of_psi_bilinear_matrix()
|
call wf_of_psi_bilinear_matrix(.True.)
|
||||||
call save_wavefunction
|
call save_wavefunction
|
||||||
|
|
||||||
deallocate (iorder, norm_sort)
|
deallocate (iorder, norm_sort)
|
||||||
|
@ -43,10 +43,12 @@ subroutine get_excitation_degree(key1,key2,degree,Nint)
|
|||||||
degree = sum(popcnt(xorvec(1:8)))
|
degree = sum(popcnt(xorvec(1:8)))
|
||||||
|
|
||||||
case default
|
case default
|
||||||
do l=1,ishft(Nint,1)
|
integer :: lmax
|
||||||
|
lmax = ishft(Nint,1)
|
||||||
|
do l=1,lmax
|
||||||
xorvec(l) = xor( key1(l), key2(l))
|
xorvec(l) = xor( key1(l), key2(l))
|
||||||
enddo
|
enddo
|
||||||
degree = sum(popcnt(xorvec(1:l)))
|
degree = sum(popcnt(xorvec(1:lmax)))
|
||||||
|
|
||||||
end select
|
end select
|
||||||
|
|
||||||
@ -245,12 +247,16 @@ subroutine get_double_excitation(det1,det2,exc,phase,Nint)
|
|||||||
if (j==k) then
|
if (j==k) then
|
||||||
nperm = nperm + popcnt(iand(det1(j,ispin), &
|
nperm = nperm + popcnt(iand(det1(j,ispin), &
|
||||||
iand( ibset(0_bit_kind,m-1)-1_bit_kind, &
|
iand( ibset(0_bit_kind,m-1)-1_bit_kind, &
|
||||||
ibclr(-1_bit_kind,n)+1_bit_kind ) ))
|
ibclr(-1_bit_kind,n)+1_bit_kind ) ))
|
||||||
|
! TODO iand( not(ishft(1_bit_kind,n+1))+1_bit_kind, &
|
||||||
|
! ishft(1_bit_kind,m)-1_bit_kind)))
|
||||||
else
|
else
|
||||||
nperm = nperm + popcnt(iand(det1(k,ispin), &
|
nperm = nperm + popcnt(iand(det1(k,ispin), &
|
||||||
ibset(0_bit_kind,m-1)-1_bit_kind))
|
ibset(0_bit_kind,m-1)-1_bit_kind))
|
||||||
|
! TODO ishft(1_bit_kind,m)-1_bit_kind))
|
||||||
if (n < bit_kind_size) then
|
if (n < bit_kind_size) then
|
||||||
nperm = nperm + popcnt(iand(det1(j,ispin), ibclr(-1_bit_kind,n) +1_bit_kind))
|
nperm = nperm + popcnt(iand(det1(j,ispin), ibclr(-1_bit_kind,n) +1_bit_kind))
|
||||||
|
! TODO ishft(1_bit_kind,m)-1_bit_kind))
|
||||||
endif
|
endif
|
||||||
do i=j+1,k-1
|
do i=j+1,k-1
|
||||||
nperm = nperm + popcnt(det1(i,ispin))
|
nperm = nperm + popcnt(det1(i,ispin))
|
||||||
@ -365,8 +371,12 @@ subroutine get_mono_excitation(det1,det2,exc,phase,Nint)
|
|||||||
if (j==k) then
|
if (j==k) then
|
||||||
nperm = popcnt(iand(det1(j,ispin), &
|
nperm = popcnt(iand(det1(j,ispin), &
|
||||||
iand(ibset(0_bit_kind,m-1)-1_bit_kind,ibclr(-1_bit_kind,n)+1_bit_kind)))
|
iand(ibset(0_bit_kind,m-1)-1_bit_kind,ibclr(-1_bit_kind,n)+1_bit_kind)))
|
||||||
|
!TODO iand( not(ishft(1_bit_kind,n+1))+1_bit_kind, &
|
||||||
|
! ishft(1_bit_kind,m)-1_bit_kind)))
|
||||||
else
|
else
|
||||||
nperm = nperm + popcnt(iand(det1(k,ispin),ibset(0_bit_kind,m-1)-1_bit_kind))
|
nperm = nperm + popcnt(iand(det1(k,ispin),ibset(0_bit_kind,m-1)-1_bit_kind))
|
||||||
|
!TODO nperm = popcnt(iand(det1(k,ispin), ishft(1_bit_kind,m)-1_bit_kind)) + &
|
||||||
|
! popcnt(iand(det1(j,ispin), not(ishft(1_bit_kind,n+1))+1_bit_kind))
|
||||||
if (n < bit_kind_size) then
|
if (n < bit_kind_size) then
|
||||||
nperm = nperm + popcnt(iand(det1(j,ispin),ibclr(-1_bit_kind,n)+1_bit_kind))
|
nperm = nperm + popcnt(iand(det1(j,ispin),ibclr(-1_bit_kind,n)+1_bit_kind))
|
||||||
endif
|
endif
|
||||||
@ -1495,14 +1505,16 @@ subroutine get_excitation_degree_vector_double_alpha_beta(key1,key2,degree,Nint,
|
|||||||
!DIR$ LOOP COUNT (1000)
|
!DIR$ LOOP COUNT (1000)
|
||||||
do i=1,sze
|
do i=1,sze
|
||||||
d = 0
|
d = 0
|
||||||
|
degree_alpha = 0
|
||||||
|
degree_beta = 0
|
||||||
!DIR$ LOOP COUNT MIN(4)
|
!DIR$ LOOP COUNT MIN(4)
|
||||||
do m=1,Nint
|
do m=1,Nint
|
||||||
d = d + popcnt(xor( key1(m,1,i), key2(m,1))) &
|
d = d + popcnt(xor( key1(m,1,i), key2(m,1))) &
|
||||||
+ popcnt(xor( key1(m,2,i), key2(m,2)))
|
+ popcnt(xor( key1(m,2,i), key2(m,2)))
|
||||||
key_tmp(m,1) = xor(key1(m,1,i),key2(m,1))
|
key_tmp(m,1) = xor(key1(m,1,i),key2(m,1))
|
||||||
key_tmp(m,2) = xor(key1(m,2,i),key2(m,2))
|
key_tmp(m,2) = xor(key1(m,2,i),key2(m,2))
|
||||||
degree_alpha = popcnt(key_tmp(m,1))
|
degree_alpha += popcnt(key_tmp(m,1))
|
||||||
degree_beta = popcnt(key_tmp(m,2))
|
degree_beta += popcnt(key_tmp(m,2))
|
||||||
enddo
|
enddo
|
||||||
if(degree_alpha .gt.3 .or. degree_beta .gt.3 )cycle !! no double excitations of same spin
|
if(degree_alpha .gt.3 .or. degree_beta .gt.3 )cycle !! no double excitations of same spin
|
||||||
degree(l) = ishft(d,-1)
|
degree(l) = ishft(d,-1)
|
||||||
@ -1653,12 +1665,13 @@ subroutine get_excitation_degree_vector_mono_or_exchange_verbose(key1,key2,degre
|
|||||||
do i=1,sze
|
do i=1,sze
|
||||||
d = 0
|
d = 0
|
||||||
exchange_1 = 0
|
exchange_1 = 0
|
||||||
|
exchange_2 = 0
|
||||||
!DIR$ LOOP COUNT MIN(4)
|
!DIR$ LOOP COUNT MIN(4)
|
||||||
do m=1,Nint
|
do m=1,Nint
|
||||||
d = d + popcnt(xor( key1(m,1,i), key2(m,1))) &
|
d = d + popcnt(xor( key1(m,1,i), key2(m,1))) &
|
||||||
+ popcnt(xor( key1(m,2,i), key2(m,2)))
|
+ popcnt(xor( key1(m,2,i), key2(m,2)))
|
||||||
exchange_1 = popcnt(xor(iand(key1(m,1,i),key1(m,2,i)),iand(key2(m,1),key2(m,2))))
|
exchange_1 += popcnt(xor(iand(key1(m,1,i),key1(m,2,i)),iand(key2(m,1),key2(m,2))))
|
||||||
exchange_2 = popcnt(iand(xor(key1(m,1,i),key2(m,1)),xor(key1(m,2,i),key2(m,2))))
|
exchange_2 += popcnt(iand(xor(key1(m,1,i),key2(m,1)),xor(key1(m,2,i),key2(m,2))))
|
||||||
enddo
|
enddo
|
||||||
if (d > 4)cycle
|
if (d > 4)cycle
|
||||||
if (d ==4)then
|
if (d ==4)then
|
||||||
@ -2217,7 +2230,7 @@ subroutine get_excitation_degree_spin(key1,key2,degree,Nint)
|
|||||||
degree = sum(popcnt(xorvec(1:4)))
|
degree = sum(popcnt(xorvec(1:4)))
|
||||||
|
|
||||||
case default
|
case default
|
||||||
do l=1,N_int
|
do l=1,Nint
|
||||||
xorvec(l) = xor( key1(l), key2(l))
|
xorvec(l) = xor( key1(l), key2(l))
|
||||||
enddo
|
enddo
|
||||||
degree = sum(popcnt(xorvec(1:Nint)))
|
degree = sum(popcnt(xorvec(1:Nint)))
|
||||||
|
@ -1217,7 +1217,6 @@ subroutine wf_of_psi_bilinear_matrix(truncate)
|
|||||||
integer :: idx
|
integer :: idx
|
||||||
integer, external :: get_index_in_psi_det_sorted_bit
|
integer, external :: get_index_in_psi_det_sorted_bit
|
||||||
double precision :: norm(N_states)
|
double precision :: norm(N_states)
|
||||||
PROVIDE psi_bilinear_matrix
|
|
||||||
|
|
||||||
do k=1,N_det
|
do k=1,N_det
|
||||||
i = psi_bilinear_matrix_rows(k)
|
i = psi_bilinear_matrix_rows(k)
|
||||||
|
12
src/Integrals_Bielec/four_idx_transform.irp.f
Normal file
12
src/Integrals_Bielec/four_idx_transform.irp.f
Normal file
@ -0,0 +1,12 @@
|
|||||||
|
program four_idx
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! 4-index transformation from AO to MO integrals
|
||||||
|
END_DOC
|
||||||
|
|
||||||
|
disk_access_mo_integrals = 'Write'
|
||||||
|
SOFT_TOUCH disk_access_mo_integrals
|
||||||
|
if (.true.) then
|
||||||
|
PROVIDE mo_bielec_integrals_in_map
|
||||||
|
endif
|
||||||
|
end
|
Loading…
Reference in New Issue
Block a user