fit of j1e in AO basis looks very different

This commit is contained in:
AbdAmmar 2024-01-17 01:59:15 +01:00
parent 430606a617
commit 3dd43d5bba
6 changed files with 182 additions and 86 deletions

View File

@ -93,7 +93,7 @@ size: (jastrow.j1e_size,nuclei.nucl_num)
type: double precision
doc: coefficients of the 1-body Jastrow in AOs
interface: ezfio
size: (nuclei.nucl_num)
size: (ao_basis.ao_num)
[j1e_expo]
type: double precision

View File

@ -1,2 +1,3 @@
nuclei
electrons
ao_basis

View File

@ -8,7 +8,11 @@ BEGIN_PROVIDER [double precision, env_val, (n_points_final_grid)]
double precision :: x, y, z, dx, dy, dz
double precision :: a, d, e, fact_r
if(env_type .eq. "Prod_Gauss") then
if(env_type .eq. "None") then
env_val = 1.d0
elseif(env_type .eq. "Prod_Gauss") then
! v(r) = \Pi_{a} [1 - \exp(-\alpha_a (r - r_a)^2)]
@ -77,7 +81,11 @@ BEGIN_PROVIDER [double precision, env_grad, (3, n_points_final_grid)]
double precision :: fact_x, fact_y, fact_z
double precision :: ax_der, ay_der, az_der, a_expo
if(env_type .eq. "Prod_Gauss") then
if(env_type .eq. "None") then
env_grad = 0.d0
elseif(env_type .eq. "Prod_Gauss") then
! v(r) = \Pi_{a} [1 - \exp(-\alpha_a (r - r_a)^2)]
@ -176,7 +184,12 @@ END_PROVIDER
PROVIDE List_env1s_square_coef List_env1s_square_expo List_env1s_square_cent
if((env_type .eq. "Prod_Gauss") .or. (env_type .eq. "Sum_Gauss")) then
if(env_type .eq. "None") then
env_square_grad = 0.d0
env_square_lapl = 0.d0
elseif((env_type .eq. "Prod_Gauss") .or. (env_type .eq. "Sum_Gauss")) then
do ipoint = 1, n_points_final_grid

View File

@ -177,29 +177,24 @@ END_PROVIDER
call get_j1e_coef_fit_ao(ao_num, coef_fit)
call ezfio_set_jastrow_j1e_coef_ao(coef_fit)
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (i, ipoint, tmp_x, tmp_y, tmp_z, &
!$OMP c) &
!$OMP SHARED (n_points_final_grid, ao_num, &
!$OMP aos_grad_in_r_array, coef_fit, &
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (i, ipoint, c) &
!$OMP SHARED (n_points_final_grid, ao_num, &
!$OMP aos_grad_in_r_array, coef_fit, &
!$OMP j1e_gradx, j1e_grady, j1e_gradz)
!$OMP DO SCHEDULE (static)
do ipoint = 1, n_points_final_grid
tmp_x = 0.d0
tmp_y = 0.d0
tmp_z = 0.d0
j1e_gradx(ipoint) = 0.d0
j1e_grady(ipoint) = 0.d0
j1e_gradz(ipoint) = 0.d0
do i = 1, ao_num
c = coef_fit(i)
tmp_x = tmp_x + c * aos_grad_in_r_array(i,ipoint,1)
tmp_y = tmp_y + c * aos_grad_in_r_array(i,ipoint,2)
tmp_z = tmp_z + c * aos_grad_in_r_array(i,ipoint,3)
j1e_gradx(ipoint) = j1e_gradx(ipoint) + c * aos_grad_in_r_array(i,ipoint,1)
j1e_grady(ipoint) = j1e_grady(ipoint) + c * aos_grad_in_r_array(i,ipoint,2)
j1e_gradz(ipoint) = j1e_gradz(ipoint) + c * aos_grad_in_r_array(i,ipoint,3)
enddo
j1e_gradx(ipoint) = tmp_x
j1e_grady(ipoint) = tmp_y
j1e_gradz(ipoint) = tmp_z
enddo
!$OMP END DO
!$OMP END PARALLEL

View File

@ -58,8 +58,8 @@ BEGIN_PROVIDER [double precision, int2_u2e_ao, (ao_num, ao_num, n_points_final_g
dy = y * env_val(ipoint)
dz = z * env_val(ipoint)
tmp0 = 0.5d0 * env_val(ipoint) * r2
tmp1 = 0.5d0 * env_val(ipoint)
tmp0 = tmp1 * r2
tmp3 = tmp_ct * env_val(ipoint)
do j = 1, ao_num
@ -124,67 +124,9 @@ BEGIN_PROVIDER [double precision, int2_grad1_u2e_ao, (ao_num, ao_num, n_points_f
if(tc_integ_type .eq. "semi-analytic") then
if((j2e_type .eq. "Mu") .and. (env_type .eq. "None")) then
PROVIDE v_ij_erf_rk_cst_mu x_v_ij_erf_rk_cst_mu
int2_grad1_u2e_ao = 0.d0
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (ipoint, i, j, x, y, z, tmp1) &
!$OMP SHARED ( ao_num, n_points_final_grid, final_grid_points &
!$OMP , v_ij_erf_rk_cst_mu, x_v_ij_erf_rk_cst_mu, int2_grad1_u2e_ao)
!$OMP DO SCHEDULE (static)
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)
do j = 1, ao_num
do i = 1, ao_num
tmp1 = v_ij_erf_rk_cst_mu(i,j,ipoint)
int2_grad1_u2e_ao(i,j,ipoint,1) = 0.5d0 * (tmp1 * x - x_v_ij_erf_rk_cst_mu(i,j,ipoint,1))
int2_grad1_u2e_ao(i,j,ipoint,2) = 0.5d0 * (tmp1 * y - x_v_ij_erf_rk_cst_mu(i,j,ipoint,2))
int2_grad1_u2e_ao(i,j,ipoint,3) = 0.5d0 * (tmp1 * z - x_v_ij_erf_rk_cst_mu(i,j,ipoint,3))
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
elseif((j2e_type .eq. "Mu") .and. (env_type .eq. "Prod_Gauss")) then
PROVIDE env_type env_val env_grad
PROVIDE v_ij_erf_rk_cst_mu_env v_ij_u_cst_mu_env_an x_v_ij_erf_rk_cst_mu_env
int2_grad1_u2e_ao = 0.d0
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (ipoint, i, j, x, y, z, tmp0, tmp1, tmp2, tmp0_x, tmp0_y, tmp0_z) &
!$OMP SHARED (ao_num, n_points_final_grid, final_grid_points, env_val, env_grad, &
!$OMP v_ij_erf_rk_cst_mu_env, v_ij_u_cst_mu_env_an, x_v_ij_erf_rk_cst_mu_env, int2_grad1_u2e_ao)
!$OMP DO SCHEDULE (static)
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)
tmp0_x = env_grad(1,ipoint)
tmp0_y = env_grad(2,ipoint)
tmp0_z = env_grad(3,ipoint)
do j = 1, ao_num
do i = 1, ao_num
tmp1 = tmp0 * v_ij_erf_rk_cst_mu_env(i,j,ipoint)
tmp2 = v_ij_u_cst_mu_env_an(i,j,ipoint)
int2_grad1_u2e_ao(i,j,ipoint,1) = tmp1 * x - tmp0 * x_v_ij_erf_rk_cst_mu_env(i,j,ipoint,1) - tmp2 * tmp0_x
int2_grad1_u2e_ao(i,j,ipoint,2) = tmp1 * y - tmp0 * x_v_ij_erf_rk_cst_mu_env(i,j,ipoint,2) - tmp2 * tmp0_y
int2_grad1_u2e_ao(i,j,ipoint,3) = tmp1 * z - tmp0 * x_v_ij_erf_rk_cst_mu_env(i,j,ipoint,3) - tmp2 * tmp0_z
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
elseif((j2e_type .eq. "Mu") .and. (env_type .eq. "Sum_Gauss")) then
if( (j2e_type .eq. "Mu") .and. &
( (env_type .eq. "None") .or. (env_type .eq. "Prod_Gauss") .or. (env_type .eq. "Sum_Gauss") ) ) then
PROVIDE mu_erf
PROVIDE env_type env_val env_grad
@ -193,8 +135,6 @@ BEGIN_PROVIDER [double precision, int2_grad1_u2e_ao, (ao_num, ao_num, n_points_f
tmp_ct = 0.5d0 / (dsqrt(dacos(-1.d0)) * mu_erf)
int2_grad1_u2e_ao = 0.d0
!$OMP PARALLEL &
!$OMP DEFAULT (NONE) &
!$OMP PRIVATE (ipoint, i, j, x, y, z, r2, dx, dy, dz, tmp1, tmp2, &
@ -300,8 +240,8 @@ subroutine get_j1e_coef_fit_ao(dim_fit, coef_fit)
allocate(u1e_tmp(n_points_final_grid))
g = 0.5d0 * (dble(elec_num) - 1.d0) / dble(elec_num)
call dgemv("T", ao_num*ao_num, n_points_final_grid, g, int2_u2e_ao(1,1,1), ao_num*ao_num, Pt, 1, 0.d0, u1e_tmp, 1)
g = -0.5d0 * (dble(elec_num) - 1.d0) / dble(elec_num)
call dgemv("T", ao_num*ao_num, n_points_final_grid, g, int2_u2e_ao, ao_num*ao_num, Pt, 1, 0.d0, u1e_tmp, 1)
FREE int2_u2e_ao
@ -340,6 +280,19 @@ subroutine get_j1e_coef_fit_ao(dim_fit, coef_fit)
! coef_fit = A_inv x b
call dgemv("N", ao_num, ao_num, 1.d0, A_inv, ao_num, b, 1, 0.d0, coef_fit, 1)
!integer :: j, k
!double precision :: tmp
!print *, ' check A_inv'
!do i = 1, ao_num
! tmp = 0.d0
! do j = 1, ao_num
! tmp += ao_overlap(i,j) * coef_fit(j)
! enddo
! tmp = tmp - b(i)
! print*, i, tmp
!enddo
deallocate(A_inv, b)
return

View File

@ -19,6 +19,12 @@ program test_non_h
touch my_extra_grid_becke my_n_pt_r_extra_grid my_n_pt_a_extra_grid
endif
PROVIDE j2e_type
PROVIDE j1e_type
PROVIDE env_type
print *, ' j2e_type = ', j2e_type
print *, ' j1e_type = ', j1e_type
print *, ' env_type = ', env_type
!call routine_fit()
@ -29,7 +35,9 @@ program test_non_h
!call test_int2_grad1_u12_square_ao()
!call test_int2_grad1_u12_ao()
call test_j1e_grad()
!call test_j1e_grad()
call test_j1e_fit_ao()
end
! ---
@ -715,3 +723,129 @@ end
! ---
subroutine test_j1e_fit_ao()
implicit none
integer :: i, j, ipoint
double precision :: g, c
double precision :: x_loops, x_dgemm, diff, thr, accu, norm
double precision, allocatable :: pa(:,:), Pb(:,:), Pt(:,:)
double precision, allocatable :: x(:), y(:), z(:)
double precision, allocatable :: x_fit(:), y_fit(:), z_fit(:), coef_fit(:)
PROVIDE mo_coef
PROVIDE int2_grad1_u2e_ao
! ---
allocate(Pa(ao_num,ao_num), Pb(ao_num,ao_num), Pt(ao_num,ao_num))
call dgemm( 'N', 'T', ao_num, ao_num, elec_alpha_num, 1.d0 &
, mo_coef, size(mo_coef, 1), mo_coef, size(mo_coef, 1) &
, 0.d0, Pa, size(Pa, 1))
if(elec_alpha_num .eq. elec_beta_num) then
Pb = Pa
else
call dgemm( 'N', 'T', ao_num, ao_num, elec_beta_num, 1.d0 &
, mo_coef, size(mo_coef, 1), mo_coef, size(mo_coef, 1) &
, 0.d0, Pb, size(Pb, 1))
endif
Pt = Pa + Pa
allocate(x(n_points_final_grid), y(n_points_final_grid), z(n_points_final_grid))
g = 0.5d0 * (dble(elec_num) - 1.d0) / dble(elec_num)
call dgemv("T", ao_num*ao_num, n_points_final_grid, g, int2_grad1_u2e_ao(1,1,1,1), ao_num*ao_num, Pt, 1, 0.d0, x, 1)
call dgemv("T", ao_num*ao_num, n_points_final_grid, g, int2_grad1_u2e_ao(1,1,1,2), ao_num*ao_num, Pt, 1, 0.d0, y, 1)
call dgemv("T", ao_num*ao_num, n_points_final_grid, g, int2_grad1_u2e_ao(1,1,1,3), ao_num*ao_num, Pt, 1, 0.d0, z, 1)
FREE int2_grad1_u2e_ao
deallocate(Pa, Pb, Pt)
! ---
allocate(x_fit(n_points_final_grid), y_fit(n_points_final_grid), z_fit(n_points_final_grid))
allocate(coef_fit(ao_num))
call get_j1e_coef_fit_ao(ao_num, coef_fit)
!print *, ' coef fit in AO:'
!print*, coef_fit
! !$OMP PARALLEL &
! !$OMP DEFAULT (NONE) &
! !$OMP PRIVATE (i, ipoint, c) &
! !$OMP SHARED (n_points_final_grid, ao_num, &
! !$OMP aos_grad_in_r_array, coef_fit, x_fit, y_fit, z_fit)
! !$OMP DO SCHEDULE (static)
do ipoint = 1, n_points_final_grid
x_fit(ipoint) = 0.d0
y_fit(ipoint) = 0.d0
z_fit(ipoint) = 0.d0
do i = 1, ao_num
c = coef_fit(i)
x_fit(ipoint) = x_fit(ipoint) + c * aos_grad_in_r_array(i,ipoint,1)
y_fit(ipoint) = y_fit(ipoint) + c * aos_grad_in_r_array(i,ipoint,2)
z_fit(ipoint) = z_fit(ipoint) + c * aos_grad_in_r_array(i,ipoint,3)
enddo
enddo
! !$OMP END DO
! !$OMP END PARALLEL
deallocate(coef_fit)
! ---
thr = 1d-10
norm = 0.d0
accu = 0.d0
do ipoint = 1, n_points_final_grid
x_loops = x (ipoint)
x_dgemm = x_fit(ipoint)
diff = dabs(x_loops - x_dgemm)
!if(diff .gt. thr) then
! print *, ' problem in j1e_gradx on:', ipoint
! print *, ' loops :', x_loops
! print *, ' dgemm :', x_dgemm
! stop
!endif
accu += diff
norm += dabs(x_loops)
x_loops = y (ipoint)
x_dgemm = y_fit(ipoint)
diff = dabs(x_loops - x_dgemm)
!if(diff .gt. thr) then
! print *, ' problem in j1e_grady on:', ipoint
! print *, ' loops :', x_loops
! print *, ' dgemm :', x_dgemm
! stop
!endif
accu += diff
norm += dabs(x_loops)
x_loops = z (ipoint)
x_dgemm = z_fit(ipoint)
diff = dabs(x_loops - x_dgemm)
!if(diff .gt. thr) then
! print *, ' problem in j1e_gradz on:', ipoint
! print *, ' loops :', x_loops
! print *, ' dgemm :', x_dgemm
! stop
!endif
accu += diff
norm += dabs(x_loops)
enddo
deallocate(x, y, z)
deallocate(x_fit, y_fit, z_fit)
print*, ' fit accuracy (%) = ', 100.d0 * accu / norm
end
! ---