10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-23 04:43:45 +01:00

added the possibility to select 3idx, 4-idx and 5idx

This commit is contained in:
eginer 2023-05-22 11:52:16 +02:00
parent 15dcd3d18f
commit 1d5ff0df66
6 changed files with 43 additions and 29 deletions

View File

@ -4,17 +4,21 @@ subroutine provide_all_three_ints_bi_ortho
! routine that provides all necessary three-electron integrals ! routine that provides all necessary three-electron integrals
END_DOC END_DOC
if(three_body_h_tc)then if(three_body_h_tc)then
if(three_e_3_idx_term)then
PROVIDE three_e_3_idx_direct_bi_ort three_e_3_idx_cycle_1_bi_ort three_e_3_idx_cycle_2_bi_ort PROVIDE three_e_3_idx_direct_bi_ort three_e_3_idx_cycle_1_bi_ort three_e_3_idx_cycle_2_bi_ort
PROVIDE three_e_3_idx_exch23_bi_ort three_e_3_idx_exch13_bi_ort three_e_3_idx_exch12_bi_ort PROVIDE three_e_3_idx_exch23_bi_ort three_e_3_idx_exch13_bi_ort three_e_3_idx_exch12_bi_ort
endif
if(three_e_4_idx_term)then
PROVIDE three_e_4_idx_direct_bi_ort three_e_4_idx_cycle_1_bi_ort three_e_4_idx_cycle_2_bi_ort PROVIDE three_e_4_idx_direct_bi_ort three_e_4_idx_cycle_1_bi_ort three_e_4_idx_cycle_2_bi_ort
PROVIDE three_e_4_idx_exch23_bi_ort three_e_4_idx_exch13_bi_ort three_e_4_idx_exch12_bi_ort PROVIDE three_e_4_idx_exch23_bi_ort three_e_4_idx_exch13_bi_ort three_e_4_idx_exch12_bi_ort
endif endif
if(.not.double_normal_ord)then if(.not.double_normal_ord.and.three_e_5_idx_term)then
PROVIDE three_e_5_idx_direct_bi_ort three_e_5_idx_cycle_1_bi_ort three_e_5_idx_cycle_2_bi_ort PROVIDE three_e_5_idx_direct_bi_ort three_e_5_idx_cycle_1_bi_ort three_e_5_idx_cycle_2_bi_ort
PROVIDE three_e_5_idx_exch23_bi_ort three_e_5_idx_exch13_bi_ort three_e_5_idx_exch12_bi_ort PROVIDE three_e_5_idx_exch23_bi_ort three_e_5_idx_exch13_bi_ort three_e_5_idx_exch12_bi_ort
else elseif (double_normal_ord .and. (.not. three_e_5_idx_term))then
PROVIDE normal_two_body_bi_orth PROVIDE normal_two_body_bi_orth
endif endif
endif
end end
subroutine diag_htilde_three_body_ints_bi_ort(Nint, key_i, hthree) subroutine diag_htilde_three_body_ints_bi_ort(Nint, key_i, hthree)

View File

@ -156,7 +156,7 @@ subroutine ac_tc_operator(iorb,ispin,key,hmono,htwoe,hthree,Nint,na,nb)
htwoe = htwoe + mo_bi_ortho_tc_two_e_jj(occ(i,other_spin),iorb) htwoe = htwoe + mo_bi_ortho_tc_two_e_jj(occ(i,other_spin),iorb)
enddo enddo
if(three_body_h_tc.and.elec_num.gt.2)then if(three_body_h_tc.and.elec_num.gt.2.and.three_e_3_idx_term)then
!!!!! 3-e part !!!!! 3-e part
!! same-spin/same-spin !! same-spin/same-spin
do j = 1, na do j = 1, na
@ -243,7 +243,7 @@ subroutine a_tc_operator(iorb,ispin,key,hmono,htwoe,hthree,Nint,na,nb)
htwoe= htwoe- mo_bi_ortho_tc_two_e_jj(occ(i,other_spin),iorb) htwoe= htwoe- mo_bi_ortho_tc_two_e_jj(occ(i,other_spin),iorb)
enddo enddo
if(three_body_h_tc.and.elec_num.gt.2)then if(three_body_h_tc.and.elec_num.gt.2.and.three_e_3_idx_term)then
!!!!! 3-e part !!!!! 3-e part
!! same-spin/same-spin !! same-spin/same-spin
do j = 1, na do j = 1, na

View File

@ -42,13 +42,13 @@ subroutine double_htilde_mu_mat_fock_bi_ortho(Nint, key_j, key_i, hmono, htwoe,
! opposite spin two-body ! opposite spin two-body
htwoe = mo_bi_ortho_tc_two_e(p2,p1,h2,h1) htwoe = mo_bi_ortho_tc_two_e(p2,p1,h2,h1)
if(three_body_h_tc.and.elec_num.gt.2)then if(three_body_h_tc.and.elec_num.gt.2)then
if(.not.double_normal_ord)then if(.not.double_normal_ord.and.three_e_5_idx_term)then
if(degree_i>degree_j)then if(degree_i>degree_j)then
call three_comp_two_e_elem(key_j,h1,h2,p1,p2,s1,s2,hthree) call three_comp_two_e_elem(key_j,h1,h2,p1,p2,s1,s2,hthree)
else else
call three_comp_two_e_elem(key_i,h1,h2,p1,p2,s1,s2,hthree) call three_comp_two_e_elem(key_i,h1,h2,p1,p2,s1,s2,hthree)
endif endif
elseif(double_normal_ord.and.elec_num.gt.2)then elseif(double_normal_ord)then
htwoe += normal_two_body_bi_orth(p2,h2,p1,h1) htwoe += normal_two_body_bi_orth(p2,h2,p1,h1)
endif endif
endif endif
@ -59,13 +59,13 @@ subroutine double_htilde_mu_mat_fock_bi_ortho(Nint, key_j, key_i, hmono, htwoe,
! exchange terms ! exchange terms
htwoe -= mo_bi_ortho_tc_two_e(p1,p2,h2,h1) htwoe -= mo_bi_ortho_tc_two_e(p1,p2,h2,h1)
if(three_body_h_tc.and.elec_num.gt.2)then if(three_body_h_tc.and.elec_num.gt.2)then
if(.not.double_normal_ord)then if(.not.double_normal_ord.and.three_e_5_idx_term)then
if(degree_i>degree_j)then if(degree_i>degree_j)then
call three_comp_two_e_elem(key_j,h1,h2,p1,p2,s1,s2,hthree) call three_comp_two_e_elem(key_j,h1,h2,p1,p2,s1,s2,hthree)
else else
call three_comp_two_e_elem(key_i,h1,h2,p1,p2,s1,s2,hthree) call three_comp_two_e_elem(key_i,h1,h2,p1,p2,s1,s2,hthree)
endif endif
elseif(double_normal_ord.and.elec_num.gt.2)then elseif(double_normal_ord)then
htwoe -= normal_two_body_bi_orth(h2,p1,h1,p2) htwoe -= normal_two_body_bi_orth(h2,p1,h1,p2)
htwoe += normal_two_body_bi_orth(h1,p1,h2,p2) htwoe += normal_two_body_bi_orth(h1,p1,h2,p2)
endif endif

View File

@ -106,7 +106,7 @@ subroutine get_single_excitation_from_fock_tc(key_i,key_j,h,p,spin,phase,hmono,h
htwoe -= buffer_x(i) htwoe -= buffer_x(i)
enddo enddo
hthree = 0.d0 hthree = 0.d0
if (three_body_h_tc.and.elec_num.gt.2)then if (three_body_h_tc.and.elec_num.gt.2.and.three_e_4_idx_term)then
call three_comp_fock_elem(key_i,h,p,spin,hthree) call three_comp_fock_elem(key_i,h,p,spin,hthree)
endif endif

View File

@ -9,33 +9,25 @@
implicit none implicit none
integer :: i, j integer :: i, j
double precision :: hmono,htwoe,hthree,htot double precision :: htot
PROVIDE N_int PROVIDE N_int
i = 1 i = 1
j = 1 j = 1
call htilde_mu_mat_bi_ortho(psi_det(1,1,j), psi_det(1,1,i), N_int, hmono, htwoe, hthree, htot) call htilde_mu_mat_opt_bi_ortho_tot(psi_det(1,1,j), psi_det(1,1,i), N_int, htot)
!$OMP PARALLEL DO SCHEDULE(GUIDED) DEFAULT(NONE) PRIVATE(i,j,hmono, htwoe, hthree, htot) & !$OMP PARALLEL DO SCHEDULE(GUIDED) DEFAULT(NONE) PRIVATE(i,j, htot) &
!$OMP SHARED (N_det, psi_det, N_int,htilde_matrix_elmt_bi_ortho) !$OMP SHARED (N_det, psi_det, N_int,htilde_matrix_elmt_bi_ortho)
do i = 1, N_det do i = 1, N_det
do j = 1, N_det do j = 1, N_det
! < J | Htilde | I > ! < J | Htilde | I >
call htilde_mu_mat_bi_ortho(psi_det(1,1,j), psi_det(1,1,i), N_int, hmono, htwoe, hthree, htot) call htilde_mu_mat_opt_bi_ortho_tot(psi_det(1,1,j), psi_det(1,1,i), N_int, htot)
!print *, ' hmono = ', hmono
!print *, ' htwoe = ', htwoe
!print *, ' hthree = ', hthree
htilde_matrix_elmt_bi_ortho(j,i) = htot htilde_matrix_elmt_bi_ortho(j,i) = htot
enddo enddo
enddo enddo
!$OMP END PARALLEL DO !$OMP END PARALLEL DO
! print*,'htilde_matrix_elmt_bi_ortho = '
! do i = 1, min(100,N_det)
! write(*,'(100(F16.10,X))')htilde_matrix_elmt_bi_ortho(1:min(100,N_det),i)
! enddo
END_PROVIDER END_PROVIDER

View File

@ -16,6 +16,24 @@ doc: If |true|, three-body terms are included
interface: ezfio,provider,ocaml interface: ezfio,provider,ocaml
default: True default: True
[three_e_3_idx_term]
type: logical
doc: If |true|, the diagonal 3-idx terms of the 3-e interaction are taken
interface: ezfio,provider,ocaml
default: True
[three_e_4_idx_term]
type: logical
doc: If |true|, the off-diagonal 4-idx terms of the 3-e interaction are taken
interface: ezfio,provider,ocaml
default: True
[three_e_5_idx_term]
type: logical
doc: If |true|, the off-diagonal 5-idx terms of the 3-e interaction are taken
interface: ezfio,provider,ocaml
default: True
[pure_three_body_h_tc] [pure_three_body_h_tc]
type: logical type: logical
doc: If |true|, pure triple excitation three-body terms are included doc: If |true|, pure triple excitation three-body terms are included