BEGIN_PROVIDER [ double precision, ao_nucl_elec_integral_pseudo, (ao_num_align,ao_num) ] implicit none ! Nucleus-pseudopotential interaction END_DOC write(output_monoints,*) 'Providing the pseudopotential integrals ' ao_nucl_elec_integral_pseudo = ao_nucl_elec_integral_pseudo_local + ao_nucl_elec_integral_pseudo_non_local END_PROVIDER BEGIN_PROVIDER [ double precision, ao_nucl_elec_integral_pseudo_local, (ao_num_align,ao_num)] implicit none BEGIN_DOC ! Local component of the pseudopotential END_DOC double precision :: alpha, beta, gama, delta integer :: num_A,num_B double precision :: A_center(3),B_center(3),C_center(3) integer :: power_A(3),power_B(3) integer :: i,j,k,l,n_pt_in,m double precision :: Vloc, Vpseudo double precision :: cpu_1, cpu_2, wall_1, wall_2, wall_0 integer :: thread_num PROVIDE output_monoints ao_nucl_elec_integral_pseudo_local = 0.d0 integer klocmax integer, allocatable :: n_k(:,:) double precision, allocatable :: v_k(:,:), dz_k(:,:) call ezfio_get_pseudo_integrals_klocmax(klocmax) allocate(n_k(nucl_num,klocmax),v_k(nucl_num,klocmax), dz_k(nucl_num,klocmax)) call ezfio_get_pseudo_integrals_v_k(v_k) call ezfio_get_pseudo_integrals_n_k(n_k) call ezfio_get_pseudo_integrals_dz_k(dz_k) !! Dump array integer, allocatable :: n_k_dump(:) double precision, allocatable :: v_k_dump(:), dz_k_dump(:) allocate(n_k_dump(1:klocmax), v_k_dump(1:klocmax), dz_k_dump(1:klocmax)) call wall_time(wall_1) call cpu_time(cpu_1) !$OMP PARALLEL & !$OMP DEFAULT (NONE) & !$OMP PRIVATE (i,j,k,l,m,alpha,beta,A_center,B_center,C_center,power_A,power_B, & !$OMP num_A,num_B,Z,c,n_pt_in, & !$OMP v_k_dump,n_k_dump, dz_k_dump, & !$OMP wall_0,wall_2,thread_num, output_monoints) & !$OMP SHARED (ao_num,ao_prim_num,ao_expo_ordered_transp,ao_power,ao_nucl,nucl_coord,ao_coef_normalized_ordered_transp, & !$OMP ao_nucl_elec_integral_pseudo_local,nucl_num,nucl_charge, & !$OMP klocmax,v_k,n_k, dz_k, & !$OMP wall_1) !$OMP DO SCHEDULE (guided) do j = 1, ao_num num_A = ao_nucl(j) power_A(1:3)= ao_power(j,1:3) A_center(1:3) = nucl_coord(num_A,1:3) do i = 1, ao_num num_B = ao_nucl(i) power_B(1:3)= ao_power(i,1:3) B_center(1:3) = nucl_coord(num_B,1:3) do l=1,ao_prim_num(j) alpha = ao_expo_ordered_transp(l,j) do m=1,ao_prim_num(i) beta = ao_expo_ordered_transp(m,i) double precision :: c c = 0.d0 do k = 1, nucl_num double precision :: Z Z = nucl_charge(k) C_center(1:3) = nucl_coord(k,1:3) v_k_dump = v_k(k,1:klocmax) n_k_dump = n_k(k,1:klocmax) dz_k_dump = dz_k(k,1:klocmax) c = c + Vloc(klocmax, v_k_dump,n_k_dump, dz_k_dump, & A_center,power_A,alpha,B_center,power_B,beta,C_center) enddo ao_nucl_elec_integral_pseudo_local(i,j) = ao_nucl_elec_integral_pseudo_local(i,j) + & ao_coef_normalized_ordered_transp(l,j)*ao_coef_normalized_ordered_transp(m,i)*c enddo enddo enddo call wall_time(wall_2) if (thread_num == 0) then if (wall_2 - wall_0 > 1.d0) then wall_0 = wall_2 write(output_monoints,*) 100.*float(j)/float(ao_num), '% in ', & wall_2-wall_1, 's' endif endif enddo !$OMP END DO !$OMP END PARALLEL deallocate(n_k,v_k, dz_k) deallocate(n_k_dump,v_k_dump, dz_k_dump) END_PROVIDER BEGIN_PROVIDER [ double precision, ao_nucl_elec_integral_pseudo_non_local, (ao_num_align,ao_num)] BEGIN_DOC ! interaction nuclear electron END_DOC implicit none double precision :: alpha, beta, gama, delta integer :: num_A,num_B double precision :: A_center(3),B_center(3),C_center(3) integer :: power_A(3),power_B(3) integer :: i,j,k,l,n_pt_in,m double precision :: Vloc, Vpseudo double precision :: cpu_1, cpu_2, wall_1, wall_2, wall_0 integer :: thread_num PROVIDE output_monoints ao_nucl_elec_integral_pseudo_non_local = 0.d0 integer :: kmax,lmax integer, allocatable :: n_kl(:,:,:) double precision, allocatable :: v_kl(:,:,:), dz_kl(:,:,:) call ezfio_get_pseudo_integrals_lmaxpo(lmax) call ezfio_get_pseudo_integrals_kmax(kmax) !lmax plus one -> lmax lmax = lmax - 1 allocate(n_kl(nucl_num,kmax,0:lmax), v_kl(nucl_num,kmax,0:lmax), dz_kl(nucl_num,kmax,0:lmax)) call ezfio_get_pseudo_integrals_n_kl(n_kl) call ezfio_get_pseudo_integrals_v_kl(v_kl) call ezfio_get_pseudo_integrals_dz_kl(dz_kl) !! Dump array integer, allocatable :: n_kl_dump(:,:) double precision, allocatable :: v_kl_dump(:,:), dz_kl_dump(:,:) allocate(n_kl_dump(kmax,0:lmax), v_kl_dump(kmax,0:lmax), dz_kl_dump(kmax,0:lmax)) ! _ ! / _. | _ | ! \_ (_| | (_ |_| | ! call wall_time(wall_1) call cpu_time(cpu_1) !$OMP PARALLEL & !$OMP DEFAULT (NONE) & !$OMP PRIVATE (i,j,k,l,m,alpha,beta,A_center,B_center,C_center,power_A,power_B, & !$OMP num_A,num_B,Z,c,n_pt_in, & !$OMP n_kl_dump, v_kl_dump, dz_kl_dump, & !$OMP wall_0,wall_2,thread_num, output_monoints) & !$OMP SHARED (ao_num,ao_prim_num,ao_expo_ordered_transp,ao_power,ao_nucl,nucl_coord,ao_coef_normalized_ordered_transp, & !$OMP ao_nucl_elec_integral_pseudo_non_local,nucl_num,nucl_charge, & !$OMP lmax,kmax,v_kl, dz_kl, n_kl, & !$OMP wall_1) !$OMP DO SCHEDULE (guided) do j = 1, ao_num num_A = ao_nucl(j) power_A(1:3)= ao_power(j,1:3) A_center(1:3) = nucl_coord(num_A,1:3) do i = 1, ao_num num_B = ao_nucl(i) power_B(1:3)= ao_power(i,1:3) B_center(1:3) = nucl_coord(num_B,1:3) do l=1,ao_prim_num(j) alpha = ao_expo_ordered_transp(l,j) do m=1,ao_prim_num(i) beta = ao_expo_ordered_transp(m,i) double precision :: c c = 0.d0 do k = 1, nucl_num double precision :: Z Z = nucl_charge(k) C_center(1:3) = nucl_coord(k,1:3) n_kl_dump = n_kl(k,1:kmax,0:lmax) v_kl_dump = v_kl(k,1:kmax,0:lmax) dz_kl_dump = dz_kl(k,1:kmax,0:lmax) c = c + Vpseudo(lmax,kmax,v_kl_dump,n_kl_dump,dz_kl_dump,A_center,power_A,alpha,B_center,power_B,beta,C_center) enddo ao_nucl_elec_integral_pseudo_non_local(i,j) = ao_nucl_elec_integral_pseudo_non_local(i,j) + & ao_coef_normalized_ordered_transp(l,j)*ao_coef_normalized_ordered_transp(m,i)*c enddo enddo enddo call wall_time(wall_2) if (thread_num == 0) then if (wall_2 - wall_0 > 1.d0) then wall_0 = wall_2 write(output_monoints,*) 100.*float(j)/float(ao_num), '% in ', & wall_2-wall_1, 's' endif endif enddo !$OMP END DO !$OMP END PARALLEL deallocate(n_kl,v_kl, dz_kl) deallocate(n_kl_dump,v_kl_dump, dz_kl_dump) END_PROVIDER