mirror of
https://github.com/QuantumPackage/qp2.git
synced 2024-12-21 11:03:29 +01:00
clean in PROVIDERS
This commit is contained in:
parent
727c70c0fa
commit
87b05b798b
@ -140,8 +140,6 @@ BEGIN_PROVIDER [ double precision, int2_grad1_u12_ao_transp, (ao_num, ao_num, 3,
|
|||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
FREE int2_grad1_u12_ao
|
|
||||||
|
|
||||||
endif
|
endif
|
||||||
|
|
||||||
call wall_time(wall1)
|
call wall_time(wall1)
|
||||||
@ -225,6 +223,8 @@ BEGIN_PROVIDER [ double precision, int2_grad1_u12_ao_t, (n_points_final_grid, 3,
|
|||||||
implicit none
|
implicit none
|
||||||
integer :: i, j, ipoint
|
integer :: i, j, ipoint
|
||||||
|
|
||||||
|
PROVIDE int2_grad1_u12_ao
|
||||||
|
|
||||||
do ipoint = 1, n_points_final_grid
|
do ipoint = 1, n_points_final_grid
|
||||||
do i = 1, ao_num
|
do i = 1, ao_num
|
||||||
do j = 1, ao_num
|
do j = 1, ao_num
|
||||||
|
@ -17,6 +17,8 @@ BEGIN_PROVIDER [ double precision, three_e_3_idx_direct_bi_ort, (mo_num, mo_num,
|
|||||||
integer :: i, j, m
|
integer :: i, j, m
|
||||||
double precision :: integral, wall1, wall0
|
double precision :: integral, wall1, wall0
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
three_e_3_idx_direct_bi_ort = 0.d0
|
three_e_3_idx_direct_bi_ort = 0.d0
|
||||||
print *, ' Providing the three_e_3_idx_direct_bi_ort ...'
|
print *, ' Providing the three_e_3_idx_direct_bi_ort ...'
|
||||||
call wall_time(wall0)
|
call wall_time(wall0)
|
||||||
@ -125,6 +127,8 @@ BEGIN_PROVIDER [ double precision, three_e_3_idx_cycle_2_bi_ort, (mo_num, mo_num
|
|||||||
integer :: i, j, m
|
integer :: i, j, m
|
||||||
double precision :: integral, wall1, wall0
|
double precision :: integral, wall1, wall0
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
three_e_3_idx_cycle_2_bi_ort = 0.d0
|
three_e_3_idx_cycle_2_bi_ort = 0.d0
|
||||||
print *, ' Providing the three_e_3_idx_cycle_2_bi_ort ...'
|
print *, ' Providing the three_e_3_idx_cycle_2_bi_ort ...'
|
||||||
call wall_time(wall0)
|
call wall_time(wall0)
|
||||||
@ -179,6 +183,8 @@ BEGIN_PROVIDER [ double precision, three_e_3_idx_exch23_bi_ort, (mo_num, mo_num,
|
|||||||
integer :: i, j, m
|
integer :: i, j, m
|
||||||
double precision :: integral, wall1, wall0
|
double precision :: integral, wall1, wall0
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
three_e_3_idx_exch23_bi_ort = 0.d0
|
three_e_3_idx_exch23_bi_ort = 0.d0
|
||||||
print*,'Providing the three_e_3_idx_exch23_bi_ort ...'
|
print*,'Providing the three_e_3_idx_exch23_bi_ort ...'
|
||||||
call wall_time(wall0)
|
call wall_time(wall0)
|
||||||
@ -233,6 +239,8 @@ BEGIN_PROVIDER [ double precision, three_e_3_idx_exch13_bi_ort, (mo_num, mo_num,
|
|||||||
integer :: i,j,m
|
integer :: i,j,m
|
||||||
double precision :: integral, wall1, wall0
|
double precision :: integral, wall1, wall0
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
three_e_3_idx_exch13_bi_ort = 0.d0
|
three_e_3_idx_exch13_bi_ort = 0.d0
|
||||||
print *, ' Providing the three_e_3_idx_exch13_bi_ort ...'
|
print *, ' Providing the three_e_3_idx_exch13_bi_ort ...'
|
||||||
call wall_time(wall0)
|
call wall_time(wall0)
|
||||||
@ -287,6 +295,8 @@ BEGIN_PROVIDER [ double precision, three_e_3_idx_exch12_bi_ort, (mo_num, mo_num,
|
|||||||
integer :: i, j, m
|
integer :: i, j, m
|
||||||
double precision :: integral, wall1, wall0
|
double precision :: integral, wall1, wall0
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
three_e_3_idx_exch12_bi_ort = 0.d0
|
three_e_3_idx_exch12_bi_ort = 0.d0
|
||||||
print *, ' Providing the three_e_3_idx_exch12_bi_ort ...'
|
print *, ' Providing the three_e_3_idx_exch12_bi_ort ...'
|
||||||
call wall_time(wall0)
|
call wall_time(wall0)
|
||||||
|
@ -1,21 +1,27 @@
|
|||||||
|
|
||||||
subroutine diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm,pt2_data,print_pt2)
|
! ---
|
||||||
|
|
||||||
|
subroutine diagonalize_CI_tc_bi_ortho(ndet, E_tc, norm, pt2_data, print_pt2)
|
||||||
|
|
||||||
|
BEGIN_DOC
|
||||||
|
! Replace the coefficients of the CI states by the coefficients of the
|
||||||
|
! eigenstates of the CI matrix
|
||||||
|
END_DOC
|
||||||
|
|
||||||
use selection_types
|
use selection_types
|
||||||
implicit none
|
implicit none
|
||||||
integer, intent(inout) :: ndet ! number of determinants from before
|
integer, intent(inout) :: ndet ! number of determinants from before
|
||||||
double precision, intent(inout) :: E_tc,norm ! E and norm from previous wave function
|
double precision, intent(inout) :: E_tc, norm ! E and norm from previous wave function
|
||||||
type(pt2_type) , intent(in) :: pt2_data ! PT2 from previous wave function
|
type(pt2_type) , intent(in) :: pt2_data ! PT2 from previous wave function
|
||||||
logical, intent(in) :: print_pt2
|
logical, intent(in) :: print_pt2
|
||||||
BEGIN_DOC
|
integer :: i, j
|
||||||
! Replace the coefficients of the CI states by the coefficients of the
|
double precision :: pt2_tmp, pt1_norm, rpt2_tmp, abs_pt2
|
||||||
! eigenstates of the CI matrix
|
|
||||||
END_DOC
|
pt2_tmp = pt2_data % pt2(1)
|
||||||
integer :: i,j
|
abs_pt2 = pt2_data % variance(1)
|
||||||
double precision :: pt2_tmp,pt1_norm,rpt2_tmp,abs_pt2
|
pt1_norm = pt2_data % overlap(1,1)
|
||||||
pt2_tmp = pt2_data % pt2(1)
|
rpt2_tmp = pt2_tmp/(1.d0 + pt1_norm)
|
||||||
abs_pt2 = pt2_data % variance(1)
|
|
||||||
pt1_norm = pt2_data % overlap(1,1)
|
|
||||||
rpt2_tmp = pt2_tmp/(1.d0 + pt1_norm)
|
|
||||||
print*,'*****'
|
print*,'*****'
|
||||||
print*,'New wave function information'
|
print*,'New wave function information'
|
||||||
print*,'N_det tc = ',N_det
|
print*,'N_det tc = ',N_det
|
||||||
@ -23,53 +29,61 @@ subroutine diagonalize_CI_tc_bi_ortho(ndet, E_tc,norm,pt2_data,print_pt2)
|
|||||||
print*,'eigval_right_tc = ',eigval_right_tc_bi_orth(1)
|
print*,'eigval_right_tc = ',eigval_right_tc_bi_orth(1)
|
||||||
print*,'Ndet, E_tc = ',N_det,eigval_right_tc_bi_orth(1)
|
print*,'Ndet, E_tc = ',N_det,eigval_right_tc_bi_orth(1)
|
||||||
print*,'*****'
|
print*,'*****'
|
||||||
if(print_pt2)then
|
|
||||||
print*,'*****'
|
if(print_pt2) then
|
||||||
print*,'previous wave function info'
|
print*,'*****'
|
||||||
print*,'norm(before) = ',norm
|
print*,'previous wave function info'
|
||||||
print*,'E(before) = ',E_tc
|
print*,'norm(before) = ',norm
|
||||||
print*,'PT1 norm = ',dsqrt(pt1_norm)
|
print*,'E(before) = ',E_tc
|
||||||
print*,'PT2 = ',pt2_tmp
|
print*,'PT1 norm = ',dsqrt(pt1_norm)
|
||||||
print*,'rPT2 = ',rpt2_tmp
|
print*,'PT2 = ',pt2_tmp
|
||||||
print*,'|PT2| = ',abs_pt2
|
print*,'rPT2 = ',rpt2_tmp
|
||||||
print*,'Positive PT2 = ',(pt2_tmp + abs_pt2)*0.5d0
|
print*,'|PT2| = ',abs_pt2
|
||||||
print*,'Negative PT2 = ',(pt2_tmp - abs_pt2)*0.5d0
|
print*,'Positive PT2 = ',(pt2_tmp + abs_pt2)*0.5d0
|
||||||
print*,'E(before) + PT2 = ',E_tc + pt2_tmp/norm
|
print*,'Negative PT2 = ',(pt2_tmp - abs_pt2)*0.5d0
|
||||||
print*,'E(before) +rPT2 = ',E_tc + rpt2_tmp/norm
|
print*,'E(before) + PT2 = ',E_tc + pt2_tmp/norm
|
||||||
write(*,'(A28,X,I10,X,100(F16.8,X))')'Ndet,E,E+PT2,E+RPT2,|PT2|=',ndet,E_tc ,E_tc + pt2_tmp/norm,E_tc + rpt2_tmp/norm,abs_pt2
|
print*,'E(before) +rPT2 = ',E_tc + rpt2_tmp/norm
|
||||||
print*,'*****'
|
write(*,'(A28,X,I10,X,100(F16.8,X))')'Ndet,E,E+PT2,E+RPT2,|PT2|=',ndet,E_tc ,E_tc + pt2_tmp/norm,E_tc + rpt2_tmp/norm,abs_pt2
|
||||||
|
print*,'*****'
|
||||||
endif
|
endif
|
||||||
|
|
||||||
psi_energy(1:N_states) = eigval_right_tc_bi_orth(1:N_states) - nuclear_repulsion
|
psi_energy(1:N_states) = eigval_right_tc_bi_orth(1:N_states) - nuclear_repulsion
|
||||||
psi_s2(1:N_states) = s2_eigvec_tc_bi_orth(1:N_states)
|
psi_s2(1:N_states) = s2_eigvec_tc_bi_orth(1:N_states)
|
||||||
|
|
||||||
E_tc = eigval_right_tc_bi_orth(1)
|
E_tc = eigval_right_tc_bi_orth(1)
|
||||||
norm = norm_ground_left_right_bi_orth
|
norm = norm_ground_left_right_bi_orth
|
||||||
ndet = N_det
|
ndet = N_det
|
||||||
do j=1,N_states
|
do j = 1, N_states
|
||||||
do i=1,N_det
|
do i = 1, N_det
|
||||||
psi_l_coef_bi_ortho(i,j) = leigvec_tc_bi_orth(i,j)
|
psi_l_coef_bi_ortho(i,j) = leigvec_tc_bi_orth(i,j)
|
||||||
psi_r_coef_bi_ortho(i,j) = reigvec_tc_bi_orth(i,j)
|
psi_r_coef_bi_ortho(i,j) = reigvec_tc_bi_orth(i,j)
|
||||||
psi_coef(i,j) = dabs(psi_l_coef_bi_ortho(i,j) * psi_r_coef_bi_ortho(i,j))
|
psi_coef(i,j) = dabs(psi_l_coef_bi_ortho(i,j) * psi_r_coef_bi_ortho(i,j))
|
||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
SOFT_TOUCH eigval_left_tc_bi_orth eigval_right_tc_bi_orth leigvec_tc_bi_orth reigvec_tc_bi_orth norm_ground_left_right_bi_orth
|
SOFT_TOUCH eigval_left_tc_bi_orth eigval_right_tc_bi_orth leigvec_tc_bi_orth reigvec_tc_bi_orth norm_ground_left_right_bi_orth
|
||||||
SOFT_TOUCH psi_l_coef_bi_ortho psi_r_coef_bi_ortho psi_coef psi_energy psi_s2
|
SOFT_TOUCH psi_l_coef_bi_ortho psi_r_coef_bi_ortho psi_coef psi_energy psi_s2
|
||||||
|
|
||||||
call save_tc_bi_ortho_wavefunction
|
call save_tc_bi_ortho_wavefunction()
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
subroutine print_CI_dressed(ndet, E_tc,norm,pt2_data,print_pt2)
|
! ---
|
||||||
|
|
||||||
|
subroutine print_CI_dressed(ndet, E_tc, norm, pt2_data, print_pt2)
|
||||||
|
|
||||||
|
BEGIN_DOC
|
||||||
|
! Replace the coefficients of the CI states by the coefficients of the
|
||||||
|
! eigenstates of the CI matrix
|
||||||
|
END_DOC
|
||||||
|
|
||||||
use selection_types
|
use selection_types
|
||||||
implicit none
|
implicit none
|
||||||
integer, intent(inout) :: ndet ! number of determinants from before
|
integer, intent(inout) :: ndet ! number of determinants from before
|
||||||
double precision, intent(inout) :: E_tc,norm ! E and norm from previous wave function
|
double precision, intent(inout) :: E_tc,norm ! E and norm from previous wave function
|
||||||
type(pt2_type) , intent(in) :: pt2_data ! PT2 from previous wave function
|
type(pt2_type) , intent(in) :: pt2_data ! PT2 from previous wave function
|
||||||
logical, intent(in) :: print_pt2
|
logical, intent(in) :: print_pt2
|
||||||
BEGIN_DOC
|
integer :: i, j
|
||||||
! Replace the coefficients of the CI states by the coefficients of the
|
|
||||||
! eigenstates of the CI matrix
|
|
||||||
END_DOC
|
|
||||||
integer :: i,j
|
|
||||||
print*,'*****'
|
print*,'*****'
|
||||||
print*,'New wave function information'
|
print*,'New wave function information'
|
||||||
print*,'N_det tc = ',N_det
|
print*,'N_det tc = ',N_det
|
||||||
@ -77,22 +91,25 @@ subroutine print_CI_dressed(ndet, E_tc,norm,pt2_data,print_pt2)
|
|||||||
print*,'eigval_right_tc = ',eigval_right_tc_bi_orth(1)
|
print*,'eigval_right_tc = ',eigval_right_tc_bi_orth(1)
|
||||||
print*,'Ndet, E_tc = ',N_det,eigval_right_tc_bi_orth(1)
|
print*,'Ndet, E_tc = ',N_det,eigval_right_tc_bi_orth(1)
|
||||||
print*,'*****'
|
print*,'*****'
|
||||||
if(print_pt2)then
|
|
||||||
print*,'*****'
|
if(print_pt2) then
|
||||||
print*,'previous wave function info'
|
print*,'*****'
|
||||||
print*,'norm(before) = ',norm
|
print*,'previous wave function info'
|
||||||
print*,'E(before) = ',E_tc
|
print*,'norm(before) = ',norm
|
||||||
print*,'PT1 norm = ',dsqrt(pt2_data % overlap(1,1))
|
print*,'E(before) = ',E_tc
|
||||||
print*,'E(before) + PT2 = ',E_tc + (pt2_data % pt2(1))/norm
|
print*,'PT1 norm = ',dsqrt(pt2_data % overlap(1,1))
|
||||||
print*,'PT2 = ',pt2_data % pt2(1)
|
print*,'E(before) + PT2 = ',E_tc + (pt2_data % pt2(1))/norm
|
||||||
print*,'Ndet, E_tc, E+PT2 = ',ndet,E_tc,E_tc + (pt2_data % pt2(1))/norm,dsqrt(pt2_data % overlap(1,1))
|
print*,'PT2 = ',pt2_data % pt2(1)
|
||||||
print*,'*****'
|
print*,'Ndet, E_tc, E+PT2 = ',ndet,E_tc,E_tc + (pt2_data % pt2(1))/norm,dsqrt(pt2_data % overlap(1,1))
|
||||||
|
print*,'*****'
|
||||||
endif
|
endif
|
||||||
|
|
||||||
E_tc = eigval_right_tc_bi_orth(1)
|
E_tc = eigval_right_tc_bi_orth(1)
|
||||||
norm = norm_ground_left_right_bi_orth
|
norm = norm_ground_left_right_bi_orth
|
||||||
ndet = N_det
|
ndet = N_det
|
||||||
do j=1,N_states
|
|
||||||
do i=1,N_det
|
do j = 1, N_states
|
||||||
|
do i = 1, N_det
|
||||||
psi_coef(i,j) = reigvec_tc_bi_orth(i,j)
|
psi_coef(i,j) = reigvec_tc_bi_orth(i,j)
|
||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
@ -100,3 +117,5 @@ subroutine print_CI_dressed(ndet, E_tc,norm,pt2_data,print_pt2)
|
|||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
|
@ -41,16 +41,20 @@ program fci
|
|||||||
my_n_pt_r_grid = 30
|
my_n_pt_r_grid = 30
|
||||||
my_n_pt_a_grid = 50
|
my_n_pt_a_grid = 50
|
||||||
touch my_grid_becke my_n_pt_r_grid my_n_pt_a_grid
|
touch my_grid_becke my_n_pt_r_grid my_n_pt_a_grid
|
||||||
|
|
||||||
pruning = -1.d0
|
pruning = -1.d0
|
||||||
touch pruning
|
touch pruning
|
||||||
|
|
||||||
! pt2_relative_error = 0.01d0
|
! pt2_relative_error = 0.01d0
|
||||||
! touch pt2_relative_error
|
! touch pt2_relative_error
|
||||||
call run_cipsi_tc
|
|
||||||
|
call run_cipsi_tc()
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
subroutine run_cipsi_tc
|
subroutine run_cipsi_tc()
|
||||||
|
|
||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
@ -59,19 +63,20 @@ subroutine run_cipsi_tc
|
|||||||
PROVIDE psi_det psi_coef mo_bi_ortho_tc_two_e mo_bi_ortho_tc_one_e
|
PROVIDE psi_det psi_coef mo_bi_ortho_tc_two_e mo_bi_ortho_tc_one_e
|
||||||
|
|
||||||
if(elec_alpha_num+elec_beta_num .ge. 3) then
|
if(elec_alpha_num+elec_beta_num .ge. 3) then
|
||||||
if(three_body_h_tc)then
|
if(three_body_h_tc) then
|
||||||
call provide_all_three_ints_bi_ortho()
|
call provide_all_three_ints_bi_ortho()
|
||||||
endif
|
endif
|
||||||
endif
|
endif
|
||||||
|
|
||||||
|
FREE int2_grad1_u12_ao
|
||||||
FREE int2_grad1_u12_bimo_transp int2_grad1_u12_ao_transp
|
FREE int2_grad1_u12_bimo_transp int2_grad1_u12_ao_transp
|
||||||
|
|
||||||
write(json_unit,json_array_open_fmt) 'fci_tc'
|
write(json_unit,json_array_open_fmt) 'fci_tc'
|
||||||
|
|
||||||
if (do_pt2) then
|
if (do_pt2) then
|
||||||
call run_stochastic_cipsi
|
call run_stochastic_cipsi()
|
||||||
else
|
else
|
||||||
call run_cipsi
|
call run_cipsi()
|
||||||
endif
|
endif
|
||||||
|
|
||||||
write(json_unit,json_dict_uopen_fmt)
|
write(json_unit,json_dict_uopen_fmt)
|
||||||
@ -83,12 +88,13 @@ subroutine run_cipsi_tc
|
|||||||
|
|
||||||
PROVIDE mo_bi_ortho_tc_one_e mo_bi_ortho_tc_two_e pt2_min_parallel_tasks
|
PROVIDE mo_bi_ortho_tc_one_e mo_bi_ortho_tc_two_e pt2_min_parallel_tasks
|
||||||
|
|
||||||
if(elec_alpha_num+elec_beta_num.ge.3)then
|
if(elec_alpha_num+elec_beta_num .ge. 3) then
|
||||||
if(three_body_h_tc)then
|
if(three_body_h_tc) then
|
||||||
call provide_all_three_ints_bi_ortho
|
call provide_all_three_ints_bi_ortho()
|
||||||
endif
|
endif
|
||||||
endif
|
endif
|
||||||
|
|
||||||
|
FREE int2_grad1_u12_ao
|
||||||
FREE int2_grad1_u12_bimo_transp int2_grad1_u12_ao_transp
|
FREE int2_grad1_u12_bimo_transp int2_grad1_u12_ao_transp
|
||||||
|
|
||||||
call run_slave_cipsi
|
call run_slave_cipsi
|
||||||
|
@ -70,7 +70,8 @@ BEGIN_PROVIDER [double precision, int2_grad1_u12_ao, (ao_num, ao_num, n_points_f
|
|||||||
|
|
||||||
elseif((j1b_type .eq. 3) .or. (j1b_type .eq. 4)) then
|
elseif((j1b_type .eq. 3) .or. (j1b_type .eq. 4)) then
|
||||||
|
|
||||||
PROVIDE v_1b_grad v_ij_erf_rk_cst_mu_j1b v_ij_u_cst_mu_j1b_an x_v_ij_erf_rk_cst_mu_j1b
|
PROVIDE v_1b_grad
|
||||||
|
PROVIDE v_ij_erf_rk_cst_mu_j1b v_ij_u_cst_mu_j1b_an x_v_ij_erf_rk_cst_mu_j1b
|
||||||
|
|
||||||
int2_grad1_u12_ao = 0.d0
|
int2_grad1_u12_ao = 0.d0
|
||||||
!$OMP PARALLEL &
|
!$OMP PARALLEL &
|
||||||
|
@ -14,81 +14,89 @@ subroutine diag_htilde_three_body_ints_bi_ort_slow(Nint, key_i, hthree)
|
|||||||
integer :: occ(Nint*bit_kind_size,2)
|
integer :: occ(Nint*bit_kind_size,2)
|
||||||
integer :: Ne(2),i,j,ii,jj,ispin,jspin,m,mm
|
integer :: Ne(2),i,j,ii,jj,ispin,jspin,m,mm
|
||||||
integer(bit_kind) :: key_i_core(Nint,2)
|
integer(bit_kind) :: key_i_core(Nint,2)
|
||||||
double precision :: direct_int, exchange_int
|
double precision :: direct_int, exchange_int, ref
|
||||||
double precision :: sym_3_e_int_from_6_idx_tensor
|
double precision, external :: sym_3_e_int_from_6_idx_tensor
|
||||||
double precision :: three_e_diag_parrallel_spin
|
double precision, external :: three_e_diag_parrallel_spin
|
||||||
|
|
||||||
if(core_tc_op)then
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
do i = 1, Nint
|
|
||||||
key_i_core(i,1) = xor(key_i(i,1),core_bitmask(i,1))
|
if(core_tc_op) then
|
||||||
key_i_core(i,2) = xor(key_i(i,2),core_bitmask(i,2))
|
do i = 1, Nint
|
||||||
enddo
|
key_i_core(i,1) = xor(key_i(i,1),core_bitmask(i,1))
|
||||||
call bitstring_to_list_ab(key_i_core,occ,Ne,Nint)
|
key_i_core(i,2) = xor(key_i(i,2),core_bitmask(i,2))
|
||||||
|
enddo
|
||||||
|
call bitstring_to_list_ab(key_i_core,occ,Ne,Nint)
|
||||||
else
|
else
|
||||||
call bitstring_to_list_ab(key_i,occ,Ne,Nint)
|
call bitstring_to_list_ab(key_i,occ,Ne,Nint)
|
||||||
endif
|
endif
|
||||||
|
|
||||||
hthree = 0.d0
|
hthree = 0.d0
|
||||||
|
|
||||||
if(Ne(1)+Ne(2).ge.3)then
|
if((Ne(1)+Ne(2)) .ge. 3) then
|
||||||
!! ! alpha/alpha/beta three-body
|
|
||||||
do i = 1, Ne(1)
|
|
||||||
ii = occ(i,1)
|
|
||||||
do j = i+1, Ne(1)
|
|
||||||
jj = occ(j,1)
|
|
||||||
do m = 1, Ne(2)
|
|
||||||
mm = occ(m,2)
|
|
||||||
! direct_int = three_body_ints_bi_ort(mm,jj,ii,mm,jj,ii) USES THE 6-IDX TENSOR
|
|
||||||
! exchange_int = three_body_ints_bi_ort(mm,jj,ii,mm,ii,jj) USES THE 6-IDX TENSOR
|
|
||||||
direct_int = three_e_3_idx_direct_bi_ort(mm,jj,ii) ! USES 3-IDX TENSOR
|
|
||||||
exchange_int = three_e_3_idx_exch12_bi_ort(mm,jj,ii) ! USES 3-IDX TENSOR
|
|
||||||
hthree += direct_int - exchange_int
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
|
|
||||||
! beta/beta/alpha three-body
|
! alpha/alpha/beta three-body
|
||||||
do i = 1, Ne(2)
|
do i = 1, Ne(1)
|
||||||
ii = occ(i,2)
|
ii = occ(i,1)
|
||||||
do j = i+1, Ne(2)
|
do j = i+1, Ne(1)
|
||||||
jj = occ(j,2)
|
jj = occ(j,1)
|
||||||
do m = 1, Ne(1)
|
do m = 1, Ne(2)
|
||||||
mm = occ(m,1)
|
mm = occ(m,2)
|
||||||
direct_int = three_e_3_idx_direct_bi_ort(mm,jj,ii)
|
!direct_int = three_body_ints_bi_ort(mm,jj,ii,mm,jj,ii) !uses the 6-idx tensor
|
||||||
exchange_int = three_e_3_idx_exch12_bi_ort(mm,jj,ii)
|
!exchange_int = three_body_ints_bi_ort(mm,jj,ii,mm,ii,jj) !uses the 6-idx tensor
|
||||||
hthree += direct_int - exchange_int
|
direct_int = three_e_3_idx_direct_bi_ort(mm,jj,ii) !uses 3-idx tensor
|
||||||
enddo
|
exchange_int = three_e_3_idx_exch12_bi_ort(mm,jj,ii) !uses 3-idx tensor
|
||||||
|
hthree += direct_int - exchange_int
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
enddo
|
enddo
|
||||||
enddo
|
|
||||||
|
|
||||||
! alpha/alpha/alpha three-body
|
! beta/beta/alpha three-body
|
||||||
do i = 1, Ne(1)
|
do i = 1, Ne(2)
|
||||||
ii = occ(i,1) ! 1
|
ii = occ(i,2)
|
||||||
do j = i+1, Ne(1)
|
do j = i+1, Ne(2)
|
||||||
jj = occ(j,1) ! 2
|
jj = occ(j,2)
|
||||||
do m = j+1, Ne(1)
|
do m = 1, Ne(1)
|
||||||
mm = occ(m,1) ! 3
|
mm = occ(m,1)
|
||||||
! ref = sym_3_e_int_from_6_idx_tensor(mm,jj,ii,mm,jj,ii) USES THE 6 IDX TENSOR
|
!direct_int = three_body_ints_bi_ort(mm,jj,ii,mm,jj,ii) !uses the 6-idx tensor
|
||||||
hthree += three_e_diag_parrallel_spin(mm,jj,ii) ! USES ONLY 3-IDX TENSORS
|
!exchange_int = three_body_ints_bi_ort(mm,jj,ii,mm,ii,jj) !uses the 6-idx tensor
|
||||||
enddo
|
direct_int = three_e_3_idx_direct_bi_ort(mm,jj,ii)
|
||||||
|
exchange_int = three_e_3_idx_exch12_bi_ort(mm,jj,ii)
|
||||||
|
hthree += direct_int - exchange_int
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
enddo
|
enddo
|
||||||
enddo
|
|
||||||
|
|
||||||
! beta/beta/beta three-body
|
! alpha/alpha/alpha three-body
|
||||||
do i = 1, Ne(2)
|
do i = 1, Ne(1)
|
||||||
ii = occ(i,2) ! 1
|
ii = occ(i,1) ! 1
|
||||||
do j = i+1, Ne(2)
|
do j = i+1, Ne(1)
|
||||||
jj = occ(j,2) ! 2
|
jj = occ(j,1) ! 2
|
||||||
do m = j+1, Ne(2)
|
do m = j+1, Ne(1)
|
||||||
mm = occ(m,2) ! 3
|
mm = occ(m,1) ! 3
|
||||||
! ref = sym_3_e_int_from_6_idx_tensor(mm,jj,ii,mm,jj,ii) USES THE 6 IDX TENSOR
|
!hthree += sym_3_e_int_from_6_idx_tensor(mm,jj,ii,mm,jj,ii) !uses the 6 idx tensor
|
||||||
hthree += three_e_diag_parrallel_spin(mm,jj,ii) ! USES ONLY 3-IDX TENSORS
|
hthree += three_e_diag_parrallel_spin(mm,jj,ii) !uses only 3-idx tensors
|
||||||
enddo
|
enddo
|
||||||
|
enddo
|
||||||
enddo
|
enddo
|
||||||
enddo
|
|
||||||
|
! beta/beta/beta three-body
|
||||||
|
do i = 1, Ne(2)
|
||||||
|
ii = occ(i,2) ! 1
|
||||||
|
do j = i+1, Ne(2)
|
||||||
|
jj = occ(j,2) ! 2
|
||||||
|
do m = j+1, Ne(2)
|
||||||
|
mm = occ(m,2) ! 3
|
||||||
|
!hthree += sym_3_e_int_from_6_idx_tensor(mm,jj,ii,mm,jj,ii) !uses the 6 idx tensor
|
||||||
|
hthree += three_e_diag_parrallel_spin(mm,jj,ii) !uses only 3-idx tensors
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
|
||||||
endif
|
endif
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
subroutine single_htilde_three_body_ints_bi_ort_slow(Nint, key_j, key_i, hthree)
|
subroutine single_htilde_three_body_ints_bi_ort_slow(Nint, key_j, key_i, hthree)
|
||||||
|
|
||||||
|
@ -1,3 +1,4 @@
|
|||||||
|
|
||||||
! ---
|
! ---
|
||||||
|
|
||||||
subroutine provide_all_three_ints_bi_ortho()
|
subroutine provide_all_three_ints_bi_ortho()
|
||||||
@ -25,9 +26,9 @@ subroutine provide_all_three_ints_bi_ortho()
|
|||||||
PROVIDE normal_two_body_bi_orth
|
PROVIDE normal_two_body_bi_orth
|
||||||
endif
|
endif
|
||||||
|
|
||||||
endif
|
endif
|
||||||
|
|
||||||
return
|
return
|
||||||
end
|
end
|
||||||
|
|
||||||
! ---
|
! ---
|
||||||
|
@ -1,32 +1,48 @@
|
|||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
BEGIN_PROVIDER [ double precision, ref_tc_energy_tot]
|
BEGIN_PROVIDER [ double precision, ref_tc_energy_tot]
|
||||||
&BEGIN_PROVIDER [ double precision, ref_tc_energy_1e]
|
&BEGIN_PROVIDER [ double precision, ref_tc_energy_1e]
|
||||||
&BEGIN_PROVIDER [ double precision, ref_tc_energy_2e]
|
&BEGIN_PROVIDER [ double precision, ref_tc_energy_2e]
|
||||||
&BEGIN_PROVIDER [ double precision, ref_tc_energy_3e]
|
&BEGIN_PROVIDER [ double precision, ref_tc_energy_3e]
|
||||||
implicit none
|
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! Various component of the TC energy for the reference "HF" Slater determinant
|
! Various component of the TC energy for the reference "HF" Slater determinant
|
||||||
END_DOC
|
END_DOC
|
||||||
double precision :: hmono, htwoe, htot, hthree
|
|
||||||
call diag_htilde_mu_mat_bi_ortho_slow(N_int,HF_bitmask , hmono, htwoe, htot)
|
implicit none
|
||||||
ref_tc_energy_1e = hmono
|
double precision :: hmono, htwoe, htot, hthree
|
||||||
ref_tc_energy_2e = htwoe
|
|
||||||
if(three_body_h_tc)then
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
call diag_htilde_three_body_ints_bi_ort_slow(N_int, HF_bitmask, hthree)
|
|
||||||
ref_tc_energy_3e = hthree
|
call diag_htilde_mu_mat_bi_ortho_slow(N_int, HF_bitmask, hmono, htwoe, htot)
|
||||||
else
|
|
||||||
ref_tc_energy_3e = 0.d0
|
ref_tc_energy_1e = hmono
|
||||||
endif
|
ref_tc_energy_2e = htwoe
|
||||||
ref_tc_energy_tot = ref_tc_energy_1e + ref_tc_energy_2e + ref_tc_energy_3e + nuclear_repulsion
|
|
||||||
END_PROVIDER
|
if(three_body_h_tc) then
|
||||||
|
call diag_htilde_three_body_ints_bi_ort_slow(N_int, HF_bitmask, hthree)
|
||||||
|
ref_tc_energy_3e = hthree
|
||||||
|
else
|
||||||
|
ref_tc_energy_3e = 0.d0
|
||||||
|
endif
|
||||||
|
|
||||||
|
ref_tc_energy_tot = ref_tc_energy_1e + ref_tc_energy_2e + ref_tc_energy_3e + nuclear_repulsion
|
||||||
|
|
||||||
|
END_PROVIDER
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
subroutine diag_htilde_mu_mat_fock_bi_ortho(Nint, det_in, hmono, htwoe, hthree, htot)
|
subroutine diag_htilde_mu_mat_fock_bi_ortho(Nint, det_in, hmono, htwoe, hthree, htot)
|
||||||
implicit none
|
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! Computes $\langle i|H|i \rangle$.
|
! Computes $\langle i|H|i \rangle$.
|
||||||
END_DOC
|
END_DOC
|
||||||
integer,intent(in) :: Nint
|
|
||||||
integer(bit_kind),intent(in) :: det_in(Nint,2)
|
implicit none
|
||||||
double precision, intent(out) :: hmono,htwoe,htot,hthree
|
integer, intent(in) :: Nint
|
||||||
|
integer(bit_kind), intent(in) :: det_in(Nint,2)
|
||||||
|
double precision, intent(out) :: hmono, htwoe, htot, hthree
|
||||||
|
|
||||||
integer(bit_kind) :: hole(Nint,2)
|
integer(bit_kind) :: hole(Nint,2)
|
||||||
integer(bit_kind) :: particle(Nint,2)
|
integer(bit_kind) :: particle(Nint,2)
|
||||||
@ -40,7 +56,6 @@ subroutine diag_htilde_mu_mat_fock_bi_ortho(Nint, det_in, hmono, htwoe, hthree,
|
|||||||
ASSERT (sum(popcnt(det_in(:,1))) == elec_alpha_num)
|
ASSERT (sum(popcnt(det_in(:,1))) == elec_alpha_num)
|
||||||
ASSERT (sum(popcnt(det_in(:,2))) == elec_beta_num)
|
ASSERT (sum(popcnt(det_in(:,2))) == elec_beta_num)
|
||||||
|
|
||||||
|
|
||||||
nexc(1) = 0
|
nexc(1) = 0
|
||||||
nexc(2) = 0
|
nexc(2) = 0
|
||||||
do i=1,Nint
|
do i=1,Nint
|
||||||
@ -55,15 +70,15 @@ subroutine diag_htilde_mu_mat_fock_bi_ortho(Nint, det_in, hmono, htwoe, hthree,
|
|||||||
enddo
|
enddo
|
||||||
|
|
||||||
if (nexc(1)+nexc(2) == 0) then
|
if (nexc(1)+nexc(2) == 0) then
|
||||||
hmono = ref_tc_energy_1e
|
hmono = ref_tc_energy_1e
|
||||||
htwoe = ref_tc_energy_2e
|
htwoe = ref_tc_energy_2e
|
||||||
hthree= ref_tc_energy_3e
|
hthree = ref_tc_energy_3e
|
||||||
htot = ref_tc_energy_tot
|
htot = ref_tc_energy_tot
|
||||||
return
|
return
|
||||||
endif
|
endif
|
||||||
|
|
||||||
!call debug_det(det_in,Nint)
|
!call debug_det(det_in,Nint)
|
||||||
integer :: tmp(2)
|
integer :: tmp(2)
|
||||||
!DIR$ FORCEINLINE
|
!DIR$ FORCEINLINE
|
||||||
call bitstring_to_list_ab(particle, occ_particle, tmp, Nint)
|
call bitstring_to_list_ab(particle, occ_particle, tmp, Nint)
|
||||||
ASSERT (tmp(1) == nexc(1)) ! Number of particles alpha
|
ASSERT (tmp(1) == nexc(1)) ! Number of particles alpha
|
||||||
@ -73,27 +88,31 @@ subroutine diag_htilde_mu_mat_fock_bi_ortho(Nint, det_in, hmono, htwoe, hthree,
|
|||||||
ASSERT (tmp(1) == nexc(1)) ! Number of holes alpha
|
ASSERT (tmp(1) == nexc(1)) ! Number of holes alpha
|
||||||
ASSERT (tmp(2) == nexc(2)) ! Number of holes beta
|
ASSERT (tmp(2) == nexc(2)) ! Number of holes beta
|
||||||
|
|
||||||
|
hmono = ref_tc_energy_1e
|
||||||
|
htwoe = ref_tc_energy_2e
|
||||||
|
hthree = ref_tc_energy_3e
|
||||||
|
|
||||||
det_tmp = ref_bitmask
|
det_tmp = ref_bitmask
|
||||||
hmono = ref_tc_energy_1e
|
|
||||||
htwoe = ref_tc_energy_2e
|
do ispin = 1, 2
|
||||||
hthree= ref_tc_energy_3e
|
|
||||||
do ispin=1,2
|
|
||||||
na = elec_num_tab(ispin)
|
na = elec_num_tab(ispin)
|
||||||
nb = elec_num_tab(iand(ispin,1)+1)
|
nb = elec_num_tab(iand(ispin,1)+1)
|
||||||
do i=1,nexc(ispin)
|
do i = 1, nexc(ispin)
|
||||||
!DIR$ FORCEINLINE
|
!DIR$ FORCEINLINE
|
||||||
call ac_tc_operator( occ_particle(i,ispin), ispin, det_tmp, hmono,htwoe,hthree, Nint,na,nb)
|
call ac_tc_operator(occ_particle(i,ispin), ispin, det_tmp, hmono, htwoe, hthree, Nint, na, nb)
|
||||||
!DIR$ FORCEINLINE
|
!DIR$ FORCEINLINE
|
||||||
call a_tc_operator ( occ_hole (i,ispin), ispin, det_tmp, hmono,htwoe,hthree, Nint,na,nb)
|
call a_tc_operator (occ_hole (i,ispin), ispin, det_tmp, hmono, htwoe, hthree, Nint, na, nb)
|
||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
htot = hmono+htwoe+hthree+nuclear_repulsion
|
|
||||||
|
htot = hmono + htwoe + hthree + nuclear_repulsion
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
subroutine ac_tc_operator(iorb,ispin,key,hmono,htwoe,hthree,Nint,na,nb)
|
subroutine ac_tc_operator(iorb,ispin,key,hmono,htwoe,hthree,Nint,na,nb)
|
||||||
use bitmasks
|
|
||||||
implicit none
|
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! Routine that computes one- and two-body energy corresponding
|
! Routine that computes one- and two-body energy corresponding
|
||||||
!
|
!
|
||||||
@ -105,6 +124,9 @@ subroutine ac_tc_operator(iorb,ispin,key,hmono,htwoe,hthree,Nint,na,nb)
|
|||||||
!
|
!
|
||||||
! and the quantities hmono,htwoe,hthree are INCREMENTED
|
! and the quantities hmono,htwoe,hthree are INCREMENTED
|
||||||
END_DOC
|
END_DOC
|
||||||
|
|
||||||
|
use bitmasks
|
||||||
|
implicit none
|
||||||
integer, intent(in) :: iorb, ispin, Nint
|
integer, intent(in) :: iorb, ispin, Nint
|
||||||
integer, intent(inout) :: na, nb
|
integer, intent(inout) :: na, nb
|
||||||
integer(bit_kind), intent(inout) :: key(Nint,2)
|
integer(bit_kind), intent(inout) :: key(Nint,2)
|
||||||
@ -460,14 +482,16 @@ subroutine a_tc_operator_no_3e(iorb,ispin,key,hmono,htwoe,Nint,na,nb)
|
|||||||
hmono = hmono - mo_bi_ortho_tc_one_e(iorb,iorb)
|
hmono = hmono - mo_bi_ortho_tc_one_e(iorb,iorb)
|
||||||
|
|
||||||
! Same spin
|
! Same spin
|
||||||
do i=1,na
|
do i = 1, na
|
||||||
htwoe= htwoe- mo_bi_ortho_tc_two_e_jj_anti(occ(i,ispin),iorb)
|
htwoe = htwoe- mo_bi_ortho_tc_two_e_jj_anti(occ(i,ispin),iorb)
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
! Opposite spin
|
! Opposite spin
|
||||||
do i=1,nb
|
do i = 1, 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
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
|
@ -106,6 +106,8 @@ subroutine diag_htilde_mu_mat_bi_ortho_slow(Nint, key_i, hmono, htwoe, htot)
|
|||||||
double precision :: get_mo_two_e_integral_tc_int
|
double precision :: get_mo_two_e_integral_tc_int
|
||||||
integer(bit_kind) :: key_i_core(Nint,2)
|
integer(bit_kind) :: key_i_core(Nint,2)
|
||||||
|
|
||||||
|
PROVIDE mo_bi_ortho_tc_two_e
|
||||||
|
|
||||||
! PROVIDE mo_two_e_integrals_tc_int_in_map mo_bi_ortho_tc_two_e
|
! PROVIDE mo_two_e_integrals_tc_int_in_map mo_bi_ortho_tc_two_e
|
||||||
!
|
!
|
||||||
! PROVIDE mo_integrals_erf_map core_energy nuclear_repulsion core_bitmask
|
! PROVIDE mo_integrals_erf_map core_energy nuclear_repulsion core_bitmask
|
||||||
@ -135,15 +137,6 @@ subroutine diag_htilde_mu_mat_bi_ortho_slow(Nint, key_i, hmono, htwoe, htot)
|
|||||||
ii = occ(i,ispin)
|
ii = occ(i,ispin)
|
||||||
hmono += mo_bi_ortho_tc_one_e(ii,ii)
|
hmono += mo_bi_ortho_tc_one_e(ii,ii)
|
||||||
|
|
||||||
! if(j1b_gauss .eq. 1) then
|
|
||||||
! print*,'j1b not implemented for bi ortho TC'
|
|
||||||
! print*,'stopping ....'
|
|
||||||
! stop
|
|
||||||
! !hmono += mo_j1b_gauss_hermI (ii,ii) &
|
|
||||||
! ! + mo_j1b_gauss_hermII (ii,ii) &
|
|
||||||
! ! + mo_j1b_gauss_nonherm(ii,ii)
|
|
||||||
! endif
|
|
||||||
|
|
||||||
! if(core_tc_op)then
|
! if(core_tc_op)then
|
||||||
! print*,'core_tc_op not already taken into account for bi ortho'
|
! print*,'core_tc_op not already taken into account for bi ortho'
|
||||||
! print*,'stopping ...'
|
! print*,'stopping ...'
|
||||||
|
@ -41,14 +41,21 @@ subroutine give_all_perm_for_three_e(n,l,k,m,j,i,idx_list,phase)
|
|||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
double precision function sym_3_e_int_from_6_idx_tensor(n,l,k,m,j,i)
|
! ---
|
||||||
implicit none
|
|
||||||
BEGIN_DOC
|
double precision function sym_3_e_int_from_6_idx_tensor(n, l, k, m, j, i)
|
||||||
! returns all good combinations of permutations of integrals with the good signs
|
|
||||||
!
|
BEGIN_DOC
|
||||||
! for a given (k^dagger l^dagger n^dagger m j i) <nlk|L|mji> when all indices have the same spins
|
! returns all good combinations of permutations of integrals with the good signs
|
||||||
END_DOC
|
!
|
||||||
integer, intent(in) :: n,l,k,m,j,i
|
! for a given (k^dagger l^dagger n^dagger m j i) <nlk|L|mji> when all indices have the same spins
|
||||||
|
END_DOC
|
||||||
|
|
||||||
|
implicit none
|
||||||
|
integer, intent(in) :: n, l, k, m, j, i
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
sym_3_e_int_from_6_idx_tensor = three_body_ints_bi_ort(n,l,k,m,j,i) & ! direct
|
sym_3_e_int_from_6_idx_tensor = three_body_ints_bi_ort(n,l,k,m,j,i) & ! direct
|
||||||
+ three_body_ints_bi_ort(n,l,k,j,i,m) & ! 1st cyclic permutation
|
+ three_body_ints_bi_ort(n,l,k,j,i,m) & ! 1st cyclic permutation
|
||||||
+ three_body_ints_bi_ort(n,l,k,i,m,j) & ! 2nd cyclic permutation
|
+ three_body_ints_bi_ort(n,l,k,i,m,j) & ! 2nd cyclic permutation
|
||||||
@ -56,8 +63,11 @@ double precision function sym_3_e_int_from_6_idx_tensor(n,l,k,m,j,i)
|
|||||||
- three_body_ints_bi_ort(n,l,k,i,j,m) & ! elec 2 is kept fixed
|
- three_body_ints_bi_ort(n,l,k,i,j,m) & ! elec 2 is kept fixed
|
||||||
- three_body_ints_bi_ort(n,l,k,m,i,j) ! elec 3 is kept fixed
|
- three_body_ints_bi_ort(n,l,k,m,i,j) ! elec 3 is kept fixed
|
||||||
|
|
||||||
|
return
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
double precision function direct_sym_3_e_int(n,l,k,m,j,i)
|
double precision function direct_sym_3_e_int(n,l,k,m,j,i)
|
||||||
implicit none
|
implicit none
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
@ -83,15 +93,25 @@ double precision function direct_sym_3_e_int(n,l,k,m,j,i)
|
|||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
double precision function three_e_diag_parrallel_spin(m,j,i)
|
! ---
|
||||||
implicit none
|
|
||||||
integer, intent(in) :: i,j,m
|
double precision function three_e_diag_parrallel_spin(m, j, i)
|
||||||
|
|
||||||
|
implicit none
|
||||||
|
integer, intent(in) :: i, j, m
|
||||||
|
|
||||||
|
PROVIDE mo_l_coef mo_r_coef
|
||||||
|
|
||||||
three_e_diag_parrallel_spin = three_e_3_idx_direct_bi_ort(m,j,i) ! direct
|
three_e_diag_parrallel_spin = three_e_3_idx_direct_bi_ort(m,j,i) ! direct
|
||||||
three_e_diag_parrallel_spin += three_e_3_idx_cycle_1_bi_ort(m,j,i) + three_e_3_idx_cycle_2_bi_ort(m,j,i) & ! two cyclic permutations
|
three_e_diag_parrallel_spin += three_e_3_idx_cycle_1_bi_ort(m,j,i) + three_e_3_idx_cycle_2_bi_ort(m,j,i) & ! two cyclic permutations
|
||||||
- three_e_3_idx_exch23_bi_ort(m,j,i) - three_e_3_idx_exch13_bi_ort(m,j,i) & ! two first exchange
|
- three_e_3_idx_exch23_bi_ort (m,j,i) - three_e_3_idx_exch13_bi_ort(m,j,i) & ! two first exchange
|
||||||
- three_e_3_idx_exch12_bi_ort(m,j,i) ! last exchange
|
- three_e_3_idx_exch12_bi_ort (m,j,i) ! last exchange
|
||||||
|
|
||||||
|
return
|
||||||
end
|
end
|
||||||
|
|
||||||
|
! ---
|
||||||
|
|
||||||
double precision function three_e_single_parrallel_spin(m,j,k,i)
|
double precision function three_e_single_parrallel_spin(m,j,k,i)
|
||||||
implicit none
|
implicit none
|
||||||
integer, intent(in) :: i,k,j,m
|
integer, intent(in) :: i,k,j,m
|
||||||
|
Loading…
Reference in New Issue
Block a user