1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-07-25 20:27:35 +02:00

Merge branch 'master' of gitlab.com:scemama/qp_plugins_scemama

This commit is contained in:
Anthony Scemama 2021-02-22 16:25:16 +01:00
commit 8371b5318f
5 changed files with 85 additions and 55 deletions

View File

@ -77,6 +77,7 @@ subroutine run
do k=1,ao_num
do j=l,ao_num
do i=k,ao_num
if (i==j .and. k<l) cycle
if (i>=j) then
integral = get_ao_two_e_integral(i,j,k,l,ao_integrals_map)
if (integral /= 0.d0) then

View File

@ -0,0 +1,72 @@
subroutine generate_fci_space
use bitmasks
implicit none
BEGIN_DOC
! Generates the complete FCI space
END_DOC
integer :: i, sze, ncore
integer(bit_kind) :: o(N_int,2)
integer(bit_kind) :: u, coremask
if (mo_num > 64) then
stop 'No more than 64 MOs'
endif
ncore = 0
coremask = 0_bit_kind
do i=1,mo_num
if (trim(mo_class(i)) == 'Core') then
ncore += 1
coremask = ibset(coremask,i-1)
endif
enddo
o(1,1) = iand(full_ijkl_bitmask(1),not(coremask))
o(1,2) = 0_bit_kind
call configuration_to_dets_size(o,n_det_alpha_unique,elec_alpha_num-ncore,N_int)
TOUCH n_det_alpha_unique
integer :: k,n,m, t, t1, t2
k=0
n = elec_alpha_num
m = mo_num - n
n = n
u = shiftl(1_bit_kind,n) -1
do while (u < shiftl(1_bit_kind,n+m))
if (iand(coremask, u) == coremask) then
k = k+1
psi_det_alpha_unique(1,k) = u
endif
t = ior(u,u-1)
t1 = t+1
t2 = shiftr((iand(not(t),t1)-1), trailz(u)+1)
u = ior(t1,t2)
enddo
call configuration_to_dets_size(o,n_det_beta_unique,elec_beta_num-ncore,N_int)
TOUCH n_det_beta_unique
k=0
n = elec_beta_num
m = mo_num - n
u = shiftl(1_bit_kind,n) -1
do while (u < shiftl(1_bit_kind,n+m))
if (iand(coremask, u) == coremask) then
k = k+1
psi_det_beta_unique(1,k) = u
endif
t = ior(u,u-1)
t1 = t+1
t2 = shiftr((iand(not(t),t1)-1), trailz(u)+1)
u = ior(t1,t2)
enddo
call generate_all_alpha_beta_det_products
print *, N_det
end

View File

@ -21,7 +21,11 @@ program mpn
e_pert(1) = hf_energy - e_pert(0) - nuclear_repulsion
do k=1,mp_order
! H_ij C^(k-1)
call h_s2_u_0_nstates_zmq(c_pert(1,k),s2,c_pert(1,k-1),1,N_det)
if (distributed_davidson) then
call H_S2_u_0_nstates_zmq (c_pert(1,k),s2,c_pert(1,k-1),1,N_det)
else
call H_S2_u_0_nstates_openmp(c_pert(1,k),s2,c_pert(1,k-1),1,N_det)
endif
if (k>1) then
e_pert(k) += c_pert(1,k)
endif
@ -46,57 +50,4 @@ program mpn
end
subroutine generate_fci_space
use bitmasks
implicit none
integer :: i, sze
integer(bit_kind) :: o(N_int,2)
if (mo_num > 64) then
stop 'No more than 64 MOs'
endif
o(:,1) = full_ijkl_bitmask(:)
o(:,2) = 0_bit_kind
call configuration_to_dets_size(o,n_det_alpha_unique,elec_alpha_num,N_int)
TOUCH n_det_alpha_unique
integer :: k,n,m, t, t1, t2
integer(bit_kind) :: u
k=0
n = elec_alpha_num
m = mo_num - n
u = shiftl(1_bit_kind,n) -1
do while (u < shiftl(1_bit_kind,n+m))
k = k+1
psi_det_alpha_unique(1,k) = u
t = ior(u,u-1)
t1 = t+1
t2 = shiftr((iand(not(t),t1)-1), trailz(u)+1)
u = ior(t1,t2)
enddo
call configuration_to_dets_size(o,n_det_beta_unique,elec_beta_num,N_int)
TOUCH n_det_beta_unique
k=0
n = elec_beta_num
m = mo_num - n
u = shiftl(1_bit_kind,n) -1
do while (u < shiftl(1_bit_kind,n+m))
k = k+1
psi_det_beta_unique(1,k) = u
t = ior(u,u-1)
t1 = t+1
t2 = shiftr((iand(not(t),t1)-1), trailz(u)+1)
u = ior(t1,t2)
enddo
call generate_all_alpha_beta_det_products
print *, N_det
end

View File

@ -88,7 +88,9 @@ subroutine run
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l)
!$OMP DO
do l=1,r
Yt(:,l) = 0.d0
do k=1,n
Yt(k,l) = 0.d0
enddo
do k=1,N_det
i = psi_bilinear_matrix_rows(k)
j = psi_bilinear_matrix_columns(k)

View File

@ -0,0 +1,4 @@
program hcore_guess_prog
implicit none
call hcore_guess
end