10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-20 04:02:16 +02:00

Merge branch 'garniron-microlist'

This commit is contained in:
Anthony Scemama 2016-01-05 01:38:45 +01:00
commit 387b891bfc
6 changed files with 237 additions and 36 deletions

View File

@ -10,7 +10,7 @@
#
#
[COMMON]
FC : gfortran -g -ffree-line-length-none -I . -static-libgcc
FC : gfortran -mavx -g -ffree-line-length-none -I . -static-libgcc
LAPACK_LIB : -llapack -lblas
IRPF90 : irpf90
IRPF90_FLAGS : --ninja --align=32
@ -22,7 +22,7 @@ IRPF90_FLAGS : --ninja --align=32
# 0 : Deactivate
#
[OPTION]
MODE : DEBUG ; [ OPT | PROFILE | DEBUG ] : Chooses the section below
MODE : OPT ; [ OPT | PROFILE | DEBUG ] : Chooses the section below
CACHE : 1 ; Enable cache_compile.py
OPENMP : 1 ; Append OpenMP flags

View File

@ -13,31 +13,31 @@ This file is autogenerad by
(** Keywords used to define input sections *)
type keyword =
| Ao_basis
| Determinants
| Determinants_by_hand
| Electrons
| Hartree_fock
| Integrals_bielec
| Mo_basis
| Nuclei
| Determinants
| Integrals_bielec
| Pseudo
| Perturbation
| Properties
| Hartree_fock
| Pseudo
;;
let keyword_to_string = function
| Ao_basis -> "AO basis"
| Determinants_by_hand -> "Determinants_by_hand"
| Determinants -> "Determinants"
| Electrons -> "Electrons"
| Hartree_fock -> "Hartree_fock"
| Integrals_bielec -> "Integrals_bielec"
| Mo_basis -> "MO basis"
| Nuclei -> "Molecule"
| Determinants -> "Determinants"
| Integrals_bielec -> "Integrals_bielec"
| Pseudo -> "Pseudo"
| Perturbation -> "Perturbation"
| Properties -> "Properties"
| Hartree_fock -> "Hartree_fock"
| Pseudo -> "Pseudo"
;;
@ -94,10 +94,10 @@ let get s =
f Pseudo.(read, to_rst)
| Perturbation ->
f Perturbation.(read, to_rst)
| Properties ->
f Properties.(read, to_rst)
| Hartree_fock ->
f Hartree_fock.(read, to_rst)
| Properties ->
f Properties.(read, to_rst)
end
with
| Sys_error msg -> (Printf.eprintf "Info: %s\n%!" msg ; "")
@ -139,8 +139,8 @@ let set str s =
| Integrals_bielec -> write Integrals_bielec.(of_rst, write) s
| Pseudo -> write Pseudo.(of_rst, write) s
| Perturbation -> write Perturbation.(of_rst, write) s
| Properties -> write Properties.(of_rst, write) s
| Hartree_fock -> write Hartree_fock.(of_rst, write) s
| Properties -> write Properties.(of_rst, write) s
| Electrons -> write Electrons.(of_rst, write) s
| Determinants_by_hand -> write Determinants_by_hand.(of_rst, write) s
| Nuclei -> write Nuclei.(of_rst, write) s
@ -192,8 +192,8 @@ let run check_only ezfio_filename =
Integrals_bielec ;
Pseudo ;
Perturbation ;
Properties ;
Hartree_fock ;
Properties ;
Mo_basis;
Determinants_by_hand ;
]
@ -212,7 +212,7 @@ let run check_only ezfio_filename =
match check_only with
| true -> ()
| false ->
Printf.sprintf "%s %s ; tput sgr0 2> /dev/null" editor temp_filename
Printf.sprintf "%s %s" editor temp_filename
|> Sys.command_exn
;

View File

@ -2,6 +2,8 @@ BEGIN_SHELL [ /usr/bin/env python ]
import perturbation
END_SHELL
subroutine perturb_buffer_$PERT(i_generator,buffer,buffer_size,e_2_pert_buffer,coef_pert_buffer,sum_e_2_pert,sum_norm_pert,sum_H_pert_diag,N_st,Nint,key_mask,fock_diag_tmp)
implicit none
BEGIN_DOC
@ -28,10 +30,19 @@ subroutine perturb_buffer_$PERT(i_generator,buffer,buffer_size,e_2_pert_buffer,c
integer :: N_minilist_gen
logical :: fullMatch
logical, external :: is_connected_to
integer(bit_kind), allocatable :: microlist(:,:,:), microlist_zero(:,:,:)
integer, allocatable :: idx_microlist(:), N_microlist(:), ptr_microlist(:), idx_microlist_zero(:)
integer :: mobiles(2), smallerlist
integer(bit_kind), allocatable :: microlist_gen(:,:,:)
integer, allocatable :: idx_microlist_gen(:), N_microlist_gen(:), ptr_microlist_gen(:)
allocate( minilist(Nint,2,N_det_selectors), &
minilist_gen(Nint,2,N_det_generators), &
idx_minilist(N_det_selectors) )
idx_minilist(N_det_selectors))
ASSERT (Nint > 0)
@ -40,27 +51,91 @@ subroutine perturb_buffer_$PERT(i_generator,buffer,buffer_size,e_2_pert_buffer,c
ASSERT (minval(sum_norm_pert) >= 0.d0)
ASSERT (N_st > 0)
call create_minilist(key_mask, psi_selectors, miniList, idx_miniList, N_det_selectors, N_minilist, Nint)
call create_minilist_find_previous(key_mask, psi_det_generators, miniList_gen, i_generator-1, N_minilist_gen, fullMatch, Nint)
if(fullMatch) then
deallocate( minilist, minilist_gen, idx_minilist )
return
end if
call create_minilist(key_mask, psi_selectors, minilist, idx_miniList, N_det_selectors, N_minilist, Nint)
allocate( microlist(Nint,2,N_minilist*4), &
idx_microlist(N_minilist*4), &
ptr_microlist(0:mo_tot_num*2+1), &
N_microlist(0:mo_tot_num*2) )
allocate( microlist_gen(Nint,2,N_minilist_gen*4), &
idx_microlist_gen(N_minilist_gen*4 ), &
ptr_microlist_gen(0:mo_tot_num*2+1), &
N_microlist_gen(0:mo_tot_num*2) )
if(key_mask(1,1) /= 0) then
call create_microlist(minilist, N_minilist, key_mask, microlist, idx_microlist, N_microlist, ptr_microlist, Nint)
call create_microlist(minilist_gen, N_minilist_gen, key_mask, microlist_gen, idx_microlist_gen, N_microlist_gen,ptr_microlist_gen,Nint)
allocate(microlist_zero(Nint,2,N_minilist))
allocate(idx_microlist_zero(N_minilist))
do i=0,mo_tot_num*2
do k=ptr_microlist(i),ptr_microlist(i+1)-1
idx_microlist(k) = idx_minilist(idx_microlist(k))
end do
end do
if(N_microlist(0) > 0) then
microlist_zero(:,:,1:N_microlist(0)) = microlist(:,:,1:N_microlist(0))
idx_microlist_zero(1:N_microlist(0)) = idx_microlist(1:N_microlist(0))
end if
end if
do i=1,buffer_size
if(is_connected_to(buffer(1,1,i), miniList_gen, Nint, N_minilist_gen)) then
cycle
end if
if (is_in_wavefunction(buffer(1,1,i),Nint)) then
cycle
endif
if(key_mask(1,1) /= 0) then
call getMobiles(buffer(:,:,i), key_mask, mobiles, Nint)
if(N_microlist(mobiles(1)) < N_microlist(mobiles(2))) then
smallerlist = mobiles(1)
else
smallerlist = mobiles(2)
end if
if(N_microlist_gen(smallerlist) > 0) then
if(is_connected_to(buffer(1,1,i), microlist_gen(:,:,ptr_microlist_gen(smallerlist):ptr_microlist_gen(smallerlist+1)-1), Nint, N_microlist_gen(smallerlist))) then
cycle
end if
end if
if(N_microlist_gen(0) > 0) then
if(is_connected_to(buffer(1,1,i), microlist_gen(:,:,1:ptr_microlist_gen(1)-1), Nint, N_microlist_gen(0))) then
cycle
end if
end if
if(N_microlist(smallerlist) > 0) then
microlist_zero(:,:,ptr_microlist(1):ptr_microlist(1)+N_microlist(smallerlist)-1) = microlist(:,:,ptr_microlist(smallerlist):ptr_microlist(smallerlist+1)-1)
idx_microlist_zero(ptr_microlist(1):ptr_microlist(1)+N_microlist(smallerlist)-1) = idx_microlist(ptr_microlist(smallerlist):ptr_microlist(smallerlist+1)-1)
! call merdge(microlist(:,:,:,smallerlist), idx_microlist(:,smallerlist), N_microlist(smallerlist), microlist(:,:,:,0), idx_microlist(:,0), N_microlist(0))
end if
call pt2_$PERT(psi_det_generators(1,1,i_generator),buffer(1,1,i), fock_diag_tmp, &
c_pert,e_2_pert,H_pert_diag,Nint,N_microlist(smallerlist)+N_microlist(0),n_st,microlist_zero(:,:,:),idx_microlist_zero(:),N_microlist(smallerlist)+N_microlist(0))
else
if(is_connected_to(buffer(1,1,i), miniList_gen, Nint, N_minilist_gen)) then
cycle
end if
call pt2_$PERT(psi_det_generators(1,1,i_generator),buffer(1,1,i), fock_diag_tmp, &
c_pert,e_2_pert,H_pert_diag,Nint,N_minilist,n_st,minilist,idx_minilist,N_minilist)
call pt2_$PERT(psi_det_generators(1,1,i_generator),buffer(1,1,i), fock_diag_tmp, &
c_pert,e_2_pert,H_pert_diag,Nint,N_minilist,n_st,minilist,idx_minilist,N_minilist)
end if
! call pt2_$PERT(psi_det_generators(1,1,i_generator),buffer(1,1,i), fock_diag_tmp, &
! c_pert,e_2_pert,H_pert_diag,Nint,N_minilist,n_st,minilist,idx_minilist,N_minilist)
do k = 1,N_st
e_2_pert_buffer(k,i) = e_2_pert(k)
@ -72,11 +147,11 @@ subroutine perturb_buffer_$PERT(i_generator,buffer,buffer_size,e_2_pert_buffer,c
enddo
deallocate( minilist, minilist_gen, idx_minilist )
deallocate( microlist, idx_microlist, N_microlist,ptr_microlist )
deallocate( microlist_gen, idx_microlist_gen,N_microlist_gen,ptr_microlist_gen )
end
subroutine perturb_buffer_by_mono_$PERT(i_generator,buffer,buffer_size,e_2_pert_buffer,coef_pert_buffer,sum_e_2_pert,sum_norm_pert,sum_H_pert_diag,N_st,Nint,key_mask,fock_diag_tmp)
implicit none
BEGIN_DOC

View File

@ -99,7 +99,7 @@ class H_apply(object):
deallocate(H_jj,iorder)
"""
s["size_max"] = "256"
s["size_max"] = "8192"
s["copy_buffer"] = """call copy_H_apply_buffer_to_wf
if (s2_eig) then
call make_s2_eigenfunction
@ -198,7 +198,7 @@ class H_apply(object):
!$ call omp_unset_lock(lck)
deallocate (e_2_pert_buffer, coef_pert_buffer)
"""
self.data["size_max"] = "256"
self.data["size_max"] = "8192"
self.data["initialization"] = """
PROVIDE psi_selectors_coef psi_selectors E_corr_per_selectors psi_det_sorted_bit
"""
@ -265,7 +265,7 @@ class H_apply(object):
double precision, intent(inout) :: select_max_out"""
self.data["params_post"] += ", select_max(min(i_generator,size(select_max,1)))"
self.data["size_max"] = "256"
self.data["size_max"] = "8192"
self.data["copy_buffer"] = """
call copy_H_apply_buffer_to_wf
if (s2_eig) then

View File

@ -97,25 +97,27 @@ end subroutine
subroutine $subroutine_diexcP(key_in, fs1, fh1, particl_1, fs2, fh2, particl_2, fock_diag_tmp, i_generator, iproc_in $parameters )
implicit none
integer(bit_kind), intent(in) :: key_in(N_int, 2), particl_1(N_int, 2), particl_2(N_int, 2)
double precision, intent(in) :: fock_diag_tmp(2,mo_tot_num+1)
integer(bit_kind) :: p1_mask(N_int, 2), p2_mask(N_int, 2), key_mask(N_int, 2)
integer,intent(in) :: fh1,fh2,fs1,fs2,i_generator,iproc_in
integer,intent(in) :: fs1,fs2,i_generator,iproc_in, fh1,fh2
integer(bit_kind) :: miniList(N_int, 2, N_det)
integer :: n_minilist, n_alpha, n_beta, deg(2), i, ni
$declarations
integer(bit_kind), parameter :: one = 1_bit_kind
p1_mask(:,:) = 0_bit_kind
p2_mask(:,:) = 0_bit_kind
p1_mask(ishft(fh1,-bit_kind_shift) + 1, fs1) = ishft(1,iand(fh1-1,bit_kind_size-1))
p2_mask(ishft(fh2,-bit_kind_shift) + 1, fs2) = ishft(1,iand(fh2-1,bit_kind_size-1))
p1_mask(ishft(fh1-1,-bit_kind_shift) + 1, fs1) = ishft(one,iand(fh1-1,bit_kind_size-1))
p2_mask(ishft(fh2-1,-bit_kind_shift) + 1, fs2) = ishft(one,iand(fh2-1,bit_kind_size-1))
key_mask(:,:) = key_in(:,:)
key_mask(ishft(fh1,-bit_kind_shift) + 1, fs1) -= ishft(1,iand(fh1-1,bit_kind_size-1))
key_mask(ishft(fh2,-bit_kind_shift) + 1, fs2) -= ishft(1,iand(fh2-1,bit_kind_size-1))
key_mask(ishft(fh1-1,-bit_kind_shift) + 1, fs1) -= ishft(one,iand(fh1-1,bit_kind_size-1))
key_mask(ishft(fh2-1,-bit_kind_shift) + 1, fs2) -= ishft(one,iand(fh2-1,bit_kind_size-1))
call $subroutine_diexcOrg(key_in, key_mask, p1_mask, particl_1, p2_mask, particl_2, fock_diag_tmp, i_generator, iproc_in $parameters )
end subroutine

View File

@ -98,6 +98,130 @@ subroutine filter_connected(key1,key2,Nint,sze,idx)
end
subroutine getMobiles(key,key_mask, mobiles,Nint)
use bitmasks
integer(bit_kind),intent(in) :: key(Nint,2), key_mask(Nint,2)
integer,intent(out) :: mobiles(2)
integer,intent(in) :: Nint
integer(bit_kind) :: mobileMask(Nint,2)
integer :: list(Nint*bit_kind_size), nel
do j=1,Nint
mobileMask(j,1) = xor(key(j,1), key_mask(j,1))
mobileMask(j,2) = xor(key(j,2), key_mask(j,2))
end do
call bitstring_to_list(mobileMask(:,1), list(:), nel, Nint)
if(nel == 2) then
mobiles(1) = list(1)
mobiles(2) = list(2)
else if(nel == 1) then
mobiles(1) = list(1)
call bitstring_to_list(mobileMask(:,2), list(:), nel, Nint)
mobiles(2) = list(1) + mo_tot_num
else
call bitstring_to_list(mobileMask(:,2), list(:), nel, Nint)
mobiles(1) = list(1) + mo_tot_num
mobiles(2) = list(2) + mo_tot_num
end if
end subroutine
subroutine create_microlist(minilist, N_minilist, key_mask, microlist, idx_microlist, N_microlist, ptr_microlist, Nint)
use bitmasks
integer, intent(in) :: Nint, N_minilist
integer(bit_kind), intent(in) :: minilist(Nint,2,N_minilist), key_mask(Nint,2)
integer, intent(out) :: N_microlist(0:mo_tot_num*2), ptr_microlist(0:mo_tot_num*2+1), idx_microlist(N_minilist*4)
integer(bit_kind), intent(out) :: microlist(Nint,2,N_minilist*4)
integer :: i,j,k,nt,n_element(2)
integer :: list(Nint*bit_kind_size,2), cur_microlist(0:mo_tot_num*2+1)
integer(bit_kind) :: key_mask_neg(Nint,2), mobileMask(Nint,2)
do i=1,Nint
key_mask_neg(i,1) = not(key_mask(i,1))
key_mask_neg(i,2) = not(key_mask(i,2))
end do
N_microlist(:) = 0
do i=1, N_minilist
do j=1,Nint
mobileMask(j,1) = iand(key_mask_neg(j,1), minilist(j,1,i))
mobileMask(j,2) = iand(key_mask_neg(j,2), minilist(j,2,i))
end do
call bitstring_to_list(mobileMask(:,1), list(:,1), n_element(1), Nint)
call bitstring_to_list(mobileMask(:,2), list(:,2), n_element(2), Nint)
if(n_element(1) + n_element(2) /= 4) then
N_microlist(0) = N_microlist(0) + 1
else
do j=1,n_element(1)
nt = list(j,1)
N_microlist(nt) = N_microlist(nt) + 1
end do
do j=1,n_element(2)
nt = list(j,2) + mo_tot_num
N_microlist(nt) = N_microlist(nt) + 1
end do
end if
end do
ptr_microlist(0) = 1
do i=1,mo_tot_num*2+1
ptr_microlist(i) = ptr_microlist(i-1) + N_microlist(i-1)
end do
cur_microlist(:) = ptr_microlist(:)
do i=1, N_minilist
do j=1,Nint
mobileMask(j,1) = iand(key_mask_neg(j,1), minilist(j,1,i))
mobileMask(j,2) = iand(key_mask_neg(j,2), minilist(j,2,i))
end do
call bitstring_to_list(mobileMask(:,1), list(:,1), n_element(1), Nint)
call bitstring_to_list(mobileMask(:,2), list(:,2), n_element(2), Nint)
if(n_element(1) + n_element(2) /= 4) then
idx_microlist(cur_microlist(0)) = i
microlist(:,:,cur_microlist(0)) = minilist(:,:,i)
cur_microlist(0) = cur_microlist(0) + 1
else
do j=1,n_element(1)
nt = list(j,1)
idx_microlist(cur_microlist(nt)) = i
microlist(:,:,cur_microlist(nt)) = minilist(:,:,i)
cur_microlist(nt) = cur_microlist(nt) + 1
end do
do j=1,n_element(2)
nt = list(j,2) + mo_tot_num
idx_microlist(cur_microlist(nt)) = i
microlist(:,:,cur_microlist(nt)) = minilist(:,:,i)
cur_microlist(nt) = cur_microlist(nt) + 1
end do
end if
end do
end subroutine
subroutine merdge(mic, idx_mic, N_mic, mic0, idx_mic0, N_mic0, Nint)
use bitmasks
integer(bit_kind) :: mic(Nint,2,N_mic), mic0(Nint,2,*)
integer :: idx_mic(N_mic), idx_mic0(N_mic0), N_mic, N_mic0
mic0(:,:,N_mic0+1:N_mic0+N_mic) = mic(:,:,:)
idx_mic0(N_mic0+1:N_mic0+N_mic) = idx_mic(:)
end subroutine
subroutine filter_connected_i_H_psi0(key1,key2,Nint,sze,idx)
use bitmasks
BEGIN_DOC