10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-06-20 04:02:22 +02:00
QuantumPackage/plugins/local/non_h_ints_mu/new_grad_tc_manu.irp.f
2024-01-16 19:07:20 +01:00

199 lines
7.0 KiB
Fortran

BEGIN_PROVIDER [ double precision, int2_grad1_u12_ao_test, (ao_num, ao_num, n_points_final_grid, 3)]
BEGIN_DOC
!
! !!!!!! WARNING !!!!!!!!!
!
! DEFINED WITH - SIGN
!
! FOR 3e-iontegrals this doesn't matter
!
! !!!!!! WARNING !!!!!!!!!
!
!
! int2_grad1_u12_ao_test(i,j,ipoint,:) = \int dr2 [-1 * \grad_r1 J(r1,r2)] \phi_i(r2) \phi_j(r2)
!
! where r1 = r(ipoint)
!
! if J(r1,r2) = u12:
!
! int2_grad1_u12_ao_test(i,j,ipoint,:) = 0.5 x \int dr2 [(r1 - r2) (erf(mu * r12)-1)r_12] \phi_i(r2) \phi_j(r2)
! = 0.5 * [ v_ij_erf_rk_cst_mu(i,j,ipoint) * r(:) - x_v_ij_erf_rk_cst_mu(i,j,ipoint,:) ]
!
! if J(r1,r2) = u12 x v1 x v2
!
! int2_grad1_u12_ao_test(i,j,ipoint,:) = v1 x [ 0.5 x \int dr2 [(r1 - r2) (erf(mu * r12)-1)r_12] v2 \phi_i(r2) \phi_j(r2) ]
! - \grad_1 v1 x [ \int dr2 u12 v2 \phi_i(r2) \phi_j(r2) ]
! = 0.5 env_val(ipoint) * v_ij_erf_rk_cst_mu_env(i,j,ipoint) * r(:)
! - 0.5 env_val(ipoint) * x_v_ij_erf_rk_cst_mu_env(i,j,ipoint,:)
! - env_grad[:,ipoint] * v_ij_u_cst_mu_env(i,j,ipoint)
!
!
END_DOC
implicit none
integer :: ipoint, i, j, m
double precision :: time0, time1
double precision :: x, y, z, tmp_x, tmp_y, tmp_z, tmp0, tmp1, tmp2
print*, ' providing int2_grad1_u12_ao_test ...'
call wall_time(time0)
if(read_tc_integ) then
open(unit=11, form="unformatted", file=trim(ezfio_filename)//'/work/int2_grad1_u12_ao_test', action="read")
read(11) int2_grad1_u12_ao_test
close(11)
else
if((j2e_type .eq. "Mu") .and. (env_type .eq. "Prod_Gauss")) then
do ipoint = 1, n_points_final_grid
x = final_grid_points(1,ipoint)
y = final_grid_points(2,ipoint)
z = final_grid_points(3,ipoint)
tmp0 = 0.5d0 * env_val(ipoint)
tmp_x = env_grad(1,ipoint)
tmp_y = env_grad(2,ipoint)
tmp_z = env_grad(3,ipoint)
do j = 1, ao_num
do i = 1, ao_num
tmp1 = tmp0 * v_ij_erf_rk_cst_mu_env_test(i,j,ipoint)
tmp2 = v_ij_u_cst_mu_env_test(i,j,ipoint)
int2_grad1_u12_ao_test(i,j,ipoint,1) = tmp1 * x - tmp0 * x_v_ij_erf_rk_cst_mu_env_test(i,j,ipoint,1) - tmp2 * tmp_x
int2_grad1_u12_ao_test(i,j,ipoint,2) = tmp1 * y - tmp0 * x_v_ij_erf_rk_cst_mu_env_test(i,j,ipoint,2) - tmp2 * tmp_y
int2_grad1_u12_ao_test(i,j,ipoint,3) = tmp1 * z - tmp0 * x_v_ij_erf_rk_cst_mu_env_test(i,j,ipoint,3) - tmp2 * tmp_z
enddo
enddo
enddo
else
print *, ' Error in int2_grad1_u12_ao_test: Unknown j2e_type = ', j2e_type
stop
endif ! j2e_type
endif
if(write_tc_integ.and.mpi_master) then
open(unit=11, form="unformatted", file=trim(ezfio_filename)//'/work/int2_grad1_u12_ao_test', action="write")
call ezfio_set_work_empty(.False.)
write(11) int2_grad1_u12_ao_test
close(11)
call ezfio_set_tc_keywords_io_tc_integ('Read')
endif
call wall_time(time1)
print*, ' Wall time for int2_grad1_u12_ao_test = ', time1 - time0
END_PROVIDER
! ---
BEGIN_PROVIDER [double precision, tc_grad_and_lapl_ao_test, (ao_num, ao_num, ao_num, ao_num)]
BEGIN_DOC
!
! tc_grad_and_lapl_ao_test(k,i,l,j) = < k l | -1/2 \Delta_1 u(r1,r2) - \grad_1 u(r1,r2) | ij >
!
! = 1/2 \int dr1 (phi_k(r1) \grad_r1 phi_i(r1) - phi_i(r1) \grad_r1 phi_k(r1)) . \int dr2 \grad_r1 u(r1,r2) \phi_l(r2) \phi_j(r2)
!
! This is obtained by integration by parts.
!
END_DOC
implicit none
integer :: ipoint, i, j, k, l, m
double precision :: weight1, contrib_x, contrib_y, contrib_z, tmp_x, tmp_y, tmp_z
double precision :: ao_k_r, ao_i_r, ao_i_dx, ao_i_dy, ao_i_dz
double precision :: time0, time1
double precision, allocatable :: ac_mat(:,:,:,:), b_mat(:,:,:,:)
print*, ' providing tc_grad_and_lapl_ao_test ...'
call wall_time(time0)
if(read_tc_integ) then
open(unit=11, form="unformatted", file=trim(ezfio_filename)//'/work/tc_grad_and_lapl_ao_test', action="read")
read(11) tc_grad_and_lapl_ao_test
close(11)
else
provide int2_grad1_u12_ao_test
allocate(b_mat(n_points_final_grid,ao_num,ao_num,3), ac_mat(ao_num,ao_num,ao_num,ao_num))
b_mat = 0.d0
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (i, k, ipoint, weight1, ao_i_r, ao_k_r) &
!$OMP SHARED (aos_in_r_array_transp, aos_grad_in_r_array_transp_bis, b_mat, &
!$OMP ao_num, n_points_final_grid, final_weight_at_r_vector)
!$OMP DO SCHEDULE (static)
do i = 1, ao_num
do k = 1, ao_num
do ipoint = 1, n_points_final_grid
weight1 = 0.5d0 * final_weight_at_r_vector(ipoint)
ao_i_r = aos_in_r_array_transp(ipoint,i)
ao_k_r = aos_in_r_array_transp(ipoint,k)
b_mat(ipoint,k,i,1) = weight1 * (ao_k_r * aos_grad_in_r_array_transp_bis(ipoint,i,1) - ao_i_r * aos_grad_in_r_array_transp_bis(ipoint,k,1))
b_mat(ipoint,k,i,2) = weight1 * (ao_k_r * aos_grad_in_r_array_transp_bis(ipoint,i,2) - ao_i_r * aos_grad_in_r_array_transp_bis(ipoint,k,2))
b_mat(ipoint,k,i,3) = weight1 * (ao_k_r * aos_grad_in_r_array_transp_bis(ipoint,i,3) - ao_i_r * aos_grad_in_r_array_transp_bis(ipoint,k,3))
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
ac_mat = 0.d0
do m = 1, 3
call dgemm( "N", "N", ao_num*ao_num, ao_num*ao_num, n_points_final_grid, 1.d0 &
, int2_grad1_u12_ao_test(1,1,1,m), ao_num*ao_num, b_mat(1,1,1,m), n_points_final_grid &
, 1.d0, ac_mat, ao_num*ao_num)
enddo
deallocate(b_mat)
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (i, j, k, l) &
!$OMP SHARED (ac_mat, tc_grad_and_lapl_ao_test, ao_num)
!$OMP DO SCHEDULE (static)
do j = 1, ao_num
do l = 1, ao_num
do i = 1, ao_num
do k = 1, ao_num
tc_grad_and_lapl_ao_test(k,i,l,j) = ac_mat(k,i,l,j) + ac_mat(l,j,k,i)
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
deallocate(ac_mat)
endif
if(write_tc_integ.and.mpi_master) then
open(unit=11, form="unformatted", file=trim(ezfio_filename)//'/work/tc_grad_and_lapl_ao_test', action="write")
call ezfio_set_work_empty(.False.)
write(11) tc_grad_and_lapl_ao_test
close(11)
call ezfio_set_tc_keywords_io_tc_integ('Read')
endif
call wall_time(time1)
print*, ' Wall time for tc_grad_and_lapl_ao_test (min) = ', (time1 - time0) / 60.d0
END_PROVIDER
! ---