2024-07-16 17:44:48 +02:00
|
|
|
double precision function get_ao_integ_chol(i,j,k,l)
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! CHOLESKY representation of the integral of the AO basis <ik|jl> or (ij|kl)
|
|
|
|
! i(r1) j(r1) 1/r12 k(r2) l(r2)
|
|
|
|
END_DOC
|
|
|
|
integer, intent(in) :: i,j,k,l
|
|
|
|
double precision, external :: ddot
|
|
|
|
get_ao_integ_chol = ddot(cholesky_ao_num, cholesky_ao_transp(1,i,j), 1, cholesky_ao_transp(1,k,l), 1)
|
|
|
|
|
|
|
|
end
|
|
|
|
|
2023-04-28 10:31:24 +02:00
|
|
|
BEGIN_PROVIDER [ double precision, cholesky_ao_transp, (cholesky_ao_num, ao_num, ao_num) ]
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! Transposed of the Cholesky vectors in AO basis set
|
|
|
|
END_DOC
|
|
|
|
integer :: i,j,k
|
|
|
|
do j=1,ao_num
|
|
|
|
do i=1,ao_num
|
2024-05-15 15:41:35 +02:00
|
|
|
do k=1,cholesky_ao_num
|
2023-04-28 10:31:24 +02:00
|
|
|
cholesky_ao_transp(k,i,j) = cholesky_ao(i,j,k)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
END_PROVIDER
|
|
|
|
|
2023-07-03 17:41:34 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
BEGIN_PROVIDER [ integer, cholesky_ao_num ]
|
2023-07-04 22:17:31 +02:00
|
|
|
&BEGIN_PROVIDER [ double precision, cholesky_ao, (ao_num, ao_num, 1) ]
|
2024-06-04 16:15:09 +02:00
|
|
|
use mmap_module
|
2023-07-11 17:31:58 +02:00
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! Cholesky vectors in AO basis: (ik|a):
|
|
|
|
! <ij|kl> = (ik|jl) = sum_a (ik|a).(a|jl)
|
|
|
|
!
|
|
|
|
! Last dimension of cholesky_ao is cholesky_ao_num
|
2024-05-31 20:07:29 +02:00
|
|
|
!
|
|
|
|
! https://mogp-emulator.readthedocs.io/en/latest/methods/proc/ProcPivotedCholesky.html
|
2024-06-20 17:27:21 +02:00
|
|
|
!
|
2024-05-31 20:07:29 +02:00
|
|
|
! https://doi.org/10.1016/j.apnum.2011.10.001 : Page 4, Algorithm 1
|
2024-06-20 17:27:21 +02:00
|
|
|
!
|
|
|
|
! https://www.diva-portal.org/smash/get/diva2:396223/FULLTEXT01.pdf
|
2023-07-11 17:31:58 +02:00
|
|
|
END_DOC
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-05-31 20:50:30 +02:00
|
|
|
integer*8 :: ndim8
|
|
|
|
integer :: rank
|
2024-06-02 19:03:05 +02:00
|
|
|
double precision :: tau, tau2
|
2024-06-20 13:43:46 +02:00
|
|
|
double precision, pointer :: L(:,:)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
double precision :: s
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-20 13:43:46 +02:00
|
|
|
double precision, allocatable :: D(:), Ltmp_p(:,:), Ltmp_q(:,:), D_sorted(:), Delta_col(:), Delta(:,:)
|
2024-06-02 19:03:05 +02:00
|
|
|
integer, allocatable :: addr1(:), addr2(:)
|
2024-06-17 14:57:48 +02:00
|
|
|
integer*8, allocatable :: Lset(:), Dset(:)
|
2023-07-11 22:31:51 +02:00
|
|
|
logical, allocatable :: computed(:)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-13 17:50:27 +02:00
|
|
|
integer :: i,j,k,m,p,q, dj, p2, q2, ii, jj
|
2024-06-04 16:15:09 +02:00
|
|
|
integer*8 :: i8, j8, p8, qj8, rank_max, np8
|
|
|
|
integer :: N, np, nq
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
double precision :: Dmax, Dmin, Qmax, f
|
|
|
|
double precision, external :: get_ao_two_e_integral
|
|
|
|
logical, external :: ao_two_e_integral_zero
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
double precision, external :: ao_two_e_integral
|
2024-05-31 20:50:30 +02:00
|
|
|
integer :: block_size, iblock
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-04 19:31:39 +02:00
|
|
|
double precision :: mem, mem0
|
2023-07-11 17:31:58 +02:00
|
|
|
double precision, external :: memory_of_double, memory_of_int
|
2024-05-31 20:50:30 +02:00
|
|
|
double precision, external :: memory_of_double8, memory_of_int8
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
integer, external :: getUnitAndOpen
|
2024-05-31 20:50:30 +02:00
|
|
|
integer :: iunit, ierr
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-05-31 20:50:30 +02:00
|
|
|
ndim8 = ao_num*ao_num*1_8
|
2024-06-02 19:16:56 +02:00
|
|
|
double precision :: wall0,wall1
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-04 16:15:09 +02:00
|
|
|
type(c_ptr) :: c_pointer(2)
|
|
|
|
integer :: fd(2)
|
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
PROVIDE nproc ao_cholesky_threshold do_direct_integrals qp_max_mem
|
2024-06-07 14:33:13 +02:00
|
|
|
PROVIDE nucl_coord ao_two_e_integral_schwartz
|
|
|
|
call set_multiple_levels_omp(.False.)
|
|
|
|
|
2024-06-02 19:16:56 +02:00
|
|
|
call wall_time(wall0)
|
2024-06-04 19:31:39 +02:00
|
|
|
|
|
|
|
! Will be reallocated at the end
|
2023-07-11 17:31:58 +02:00
|
|
|
deallocate(cholesky_ao)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (read_ao_cholesky) then
|
2024-06-07 16:09:53 +02:00
|
|
|
print *, 'Reading Cholesky AO vectors from disk...'
|
2023-07-11 17:31:58 +02:00
|
|
|
iunit = getUnitAndOpen(trim(ezfio_work_dir)//'cholesky_ao', 'R')
|
|
|
|
read(iunit) rank
|
|
|
|
allocate(cholesky_ao(ao_num,ao_num,rank), stat=ierr)
|
|
|
|
read(iunit) cholesky_ao
|
|
|
|
close(iunit)
|
|
|
|
cholesky_ao_num = rank
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
else
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-03-20 09:20:21 +01:00
|
|
|
call set_multiple_levels_omp(.False.)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 22:17:40 +02:00
|
|
|
if (do_direct_integrals) then
|
|
|
|
if (ao_two_e_integral(1,1,1,1) < huge(1.d0)) then
|
|
|
|
! Trigger providers inside ao_two_e_integral
|
|
|
|
continue
|
|
|
|
endif
|
|
|
|
else
|
2023-07-11 17:31:58 +02:00
|
|
|
PROVIDE ao_two_e_integrals_in_map
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
tau = ao_cholesky_threshold
|
2024-06-02 19:03:05 +02:00
|
|
|
tau2 = tau*tau
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
rank = 0
|
|
|
|
|
|
|
|
allocate( D(ndim8), Lset(ndim8), Dset(ndim8), D_sorted(ndim8))
|
2024-06-20 13:43:46 +02:00
|
|
|
allocate( addr1(ndim8), addr2(ndim8), Delta_col(ndim8), computed(ndim8) )
|
2024-06-17 14:57:48 +02:00
|
|
|
|
|
|
|
call resident_memory(mem0)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 22:17:40 +02:00
|
|
|
call print_memory_usage()
|
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
print *, ''
|
|
|
|
print *, 'Cholesky decomposition of AO integrals'
|
|
|
|
print *, '======================================'
|
|
|
|
print *, ''
|
|
|
|
print *, '============ ============='
|
2024-06-02 19:03:05 +02:00
|
|
|
print *, ' Rank Threshold'
|
2023-07-11 17:31:58 +02:00
|
|
|
print *, '============ ============='
|
2023-08-21 09:56:07 +02:00
|
|
|
|
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! 1.
|
2024-06-17 14:57:48 +02:00
|
|
|
i8=0
|
2023-07-11 17:31:58 +02:00
|
|
|
do j=1,ao_num
|
|
|
|
do i=1,ao_num
|
2024-06-17 14:57:48 +02:00
|
|
|
i8 = i8+1
|
|
|
|
addr1(i8) = i
|
|
|
|
addr2(i8) = j
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (do_direct_integrals) then
|
2024-06-17 14:57:48 +02:00
|
|
|
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(i8) SCHEDULE(dynamic,21)
|
2024-06-04 16:15:09 +02:00
|
|
|
do i8=ndim8,1,-1
|
2024-05-31 20:50:30 +02:00
|
|
|
D(i8) = ao_two_e_integral(addr1(i8), addr2(i8), &
|
|
|
|
addr1(i8), addr2(i8))
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
else
|
2024-06-17 14:57:48 +02:00
|
|
|
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(i8) SCHEDULE(dynamic,21)
|
2024-06-04 16:15:09 +02:00
|
|
|
do i8=ndim8,1,-1
|
2024-05-31 20:50:30 +02:00
|
|
|
D(i8) = get_ao_two_e_integral(addr1(i8), addr1(i8), &
|
2024-06-17 14:57:48 +02:00
|
|
|
addr2(i8), addr2(i8), ao_integrals_map)
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
endif
|
2024-06-17 14:57:48 +02:00
|
|
|
|
2024-06-04 16:15:09 +02:00
|
|
|
D_sorted(:) = -D(:)
|
|
|
|
call dsort_noidx_big(D_sorted,ndim8)
|
2024-06-17 14:57:48 +02:00
|
|
|
D_sorted(:) = -D_sorted(:)
|
2024-06-04 16:15:09 +02:00
|
|
|
Dmax = D_sorted(1)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! 2.
|
2024-06-04 19:31:39 +02:00
|
|
|
np8=0_8
|
|
|
|
do p8=1,ndim8
|
2024-06-20 13:43:46 +02:00
|
|
|
if ( Dmax*D(p8) >= tau2 ) then
|
2024-06-04 19:31:39 +02:00
|
|
|
np8 = np8+1_8
|
|
|
|
Lset(np8) = p8
|
|
|
|
endif
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2024-07-03 14:52:11 +02:00
|
|
|
if (np8 > ndim8) stop 'np>ndim8'
|
|
|
|
np = int(np8,4)
|
2024-06-17 14:57:48 +02:00
|
|
|
if (np <= 0) stop 'np<=0'
|
|
|
|
|
2024-07-16 17:44:48 +02:00
|
|
|
rank_max = np
|
2024-07-29 16:15:48 +02:00
|
|
|
! Avoid too large arrays when there are many electrons
|
|
|
|
if (elec_num > 10) then
|
|
|
|
rank_max = min(np,20*elec_num*elec_num)
|
|
|
|
endif
|
2024-06-17 14:57:48 +02:00
|
|
|
call mmap(trim(ezfio_work_dir)//'cholesky_ao_tmp', (/ ndim8, rank_max /), 8, fd(1), .False., .True., c_pointer(1))
|
|
|
|
call c_f_pointer(c_pointer(1), L, (/ ndim8, rank_max /))
|
|
|
|
|
|
|
|
! Deleting the file while it is open makes the file invisible on the filesystem,
|
|
|
|
! and automatically deleted, even if the program crashes
|
|
|
|
iunit = getUnitAndOpen(trim(ezfio_work_dir)//'cholesky_ao_tmp', 'R')
|
|
|
|
close(iunit,status='delete')
|
|
|
|
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! 3.
|
|
|
|
N = 0
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! 4.
|
|
|
|
i = 0
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
mem = memory_of_double(np) & ! Delta(np,nq)
|
|
|
|
+ (np+1)*memory_of_double(block_size) ! Ltmp_p(np,block_size) + Ltmp_q(nq,block_size)
|
|
|
|
|
|
|
|
! call check_mem(mem)
|
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! 5.
|
2024-06-17 14:57:48 +02:00
|
|
|
do while ( (Dmax > tau).and.(np > 0) )
|
2023-07-11 17:31:58 +02:00
|
|
|
! a.
|
|
|
|
i = i+1
|
2023-08-21 09:56:07 +02:00
|
|
|
|
|
|
|
|
2024-05-31 20:07:29 +02:00
|
|
|
block_size = max(N,24)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
! Determine nq so that Delta fits in memory
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
s = 0.1d0
|
|
|
|
Dmin = max(s*Dmax,tau)
|
|
|
|
do nq=2,np-1
|
|
|
|
if (D_sorted(nq) < Dmin) exit
|
|
|
|
enddo
|
2024-05-31 20:07:29 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
do while (.True.)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
mem = mem0 &
|
|
|
|
+ np*memory_of_double(nq) & ! Delta(np,nq)
|
2024-06-20 13:43:46 +02:00
|
|
|
+ (np+nq)*memory_of_double(block_size) ! Ltmp_p(np,block_size) + Ltmp_q(nq,block_size)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
if (mem > qp_max_mem*0.5d0) then
|
2024-06-20 13:43:46 +02:00
|
|
|
Dmin = D_sorted(nq/2)
|
|
|
|
do ii=nq/2,np-1
|
|
|
|
if (D_sorted(ii) < Dmin) then
|
|
|
|
nq = ii
|
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
2023-07-11 17:31:58 +02:00
|
|
|
else
|
|
|
|
exit
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2024-06-20 13:43:46 +02:00
|
|
|
!call print_memory_usage
|
|
|
|
!print *, 'np, nq, Predicted memory: ', np, nq, mem
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
if (nq <= 0) then
|
|
|
|
print *, nq
|
|
|
|
stop 'bug in cholesky: nq <= 0'
|
2023-07-11 17:31:58 +02:00
|
|
|
endif
|
2024-06-04 19:31:39 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
Dmin = D_sorted(nq)
|
|
|
|
nq=0
|
|
|
|
do p=1,np
|
|
|
|
if ( D(Lset(p)) >= Dmin ) then
|
|
|
|
nq = nq+1
|
|
|
|
Dset(nq) = Lset(p)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
allocate(Delta(np,nq))
|
2023-07-11 17:31:58 +02:00
|
|
|
allocate(Ltmp_p(np,block_size), stat=ierr)
|
2024-05-31 20:07:29 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (ierr /= 0) then
|
2023-07-11 22:17:40 +02:00
|
|
|
call print_memory_usage()
|
2023-07-11 17:31:58 +02:00
|
|
|
print *, irp_here, ': allocation failed : (Ltmp_p(np,block_size))'
|
|
|
|
stop -1
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
allocate(Ltmp_q(nq,block_size), stat=ierr)
|
2024-05-31 20:07:29 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (ierr /= 0) then
|
2023-07-11 22:17:40 +02:00
|
|
|
call print_memory_usage()
|
2023-07-11 17:31:58 +02:00
|
|
|
print *, irp_here, ': allocation failed : (Ltmp_q(nq,block_size))'
|
|
|
|
stop -1
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
|
|
|
|
2024-06-20 13:43:46 +02:00
|
|
|
computed(1:nq) = .False.
|
2023-08-21 09:56:07 +02:00
|
|
|
|
|
|
|
|
2024-06-05 02:51:12 +02:00
|
|
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(k,p,q)
|
2023-07-11 17:31:58 +02:00
|
|
|
do k=1,N
|
2024-06-02 19:03:05 +02:00
|
|
|
!$OMP DO
|
2023-07-11 17:31:58 +02:00
|
|
|
do p=1,np
|
2024-06-04 19:31:39 +02:00
|
|
|
Ltmp_p(p,k) = L(Lset(p),k)
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2024-06-02 19:03:05 +02:00
|
|
|
!$OMP END DO NOWAIT
|
|
|
|
|
|
|
|
!$OMP DO
|
2023-07-11 17:31:58 +02:00
|
|
|
do q=1,nq
|
2024-06-04 19:31:39 +02:00
|
|
|
Ltmp_q(q,k) = L(Dset(q),k)
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2024-06-02 19:03:05 +02:00
|
|
|
!$OMP END DO NOWAIT
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
|
|
|
!$OMP BARRIER
|
|
|
|
!$OMP END PARALLEL
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (N>0) then
|
2024-06-14 16:26:23 +02:00
|
|
|
|
|
|
|
call dgemm('N', 'T', np, nq, N, -1.d0, &
|
|
|
|
Ltmp_p(1,1), np, Ltmp_q(1,1), nq, 0.d0, Delta, np)
|
|
|
|
|
2024-06-04 19:31:39 +02:00
|
|
|
else
|
2024-06-14 16:26:23 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(q,j)
|
2024-06-04 19:31:39 +02:00
|
|
|
do q=1,nq
|
2024-06-17 14:57:48 +02:00
|
|
|
Delta(:,q) = 0.d0
|
2024-06-04 19:31:39 +02:00
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
2024-06-13 17:50:27 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! f.
|
|
|
|
Qmax = D(Dset(1))
|
|
|
|
do q=1,nq
|
|
|
|
Qmax = max(Qmax, D(Dset(q)))
|
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! g.
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
iblock = 0
|
|
|
|
do j=1,nq
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-05-31 20:50:30 +02:00
|
|
|
if ( (Qmax <= Dmin).or.(N+j*1_8 > ndim8) ) exit
|
2024-06-17 14:57:48 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! i.
|
|
|
|
rank = N+j
|
2024-06-17 14:57:48 +02:00
|
|
|
if (rank == rank_max) then
|
|
|
|
print *, 'cholesky: rank_max reached'
|
|
|
|
exit
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (iblock == block_size) then
|
2024-06-13 17:50:27 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
call dgemm('N','T',np,nq,block_size,-1.d0, &
|
2024-06-04 19:31:39 +02:00
|
|
|
Ltmp_p, np, Ltmp_q, nq, 1.d0, Delta, np)
|
2024-06-13 17:50:27 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
iblock = 0
|
2024-06-13 17:50:27 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! ii.
|
|
|
|
do dj=1,nq
|
2024-05-31 20:50:30 +02:00
|
|
|
qj8 = Dset(dj)
|
|
|
|
if (D(qj8) == Qmax) then
|
2023-07-11 17:31:58 +02:00
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-02 19:03:05 +02:00
|
|
|
do i8=1,ndim8
|
|
|
|
L(i8, rank) = 0.d0
|
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-04 19:31:39 +02:00
|
|
|
iblock = iblock+1
|
|
|
|
!$OMP PARALLEL DO PRIVATE(p)
|
|
|
|
do p=1,np
|
|
|
|
Ltmp_p(p,iblock) = Delta(p,dj)
|
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
|
2023-07-11 22:31:51 +02:00
|
|
|
if (.not.computed(dj)) then
|
|
|
|
m = dj
|
2024-06-02 19:03:05 +02:00
|
|
|
if (do_direct_integrals) then
|
2024-06-04 19:31:39 +02:00
|
|
|
!$OMP PARALLEL DO PRIVATE(k) SCHEDULE(dynamic,21)
|
|
|
|
do k=1,np
|
2024-06-17 14:57:48 +02:00
|
|
|
Delta_col(k) = 0.d0
|
2024-06-02 19:03:05 +02:00
|
|
|
if (.not.ao_two_e_integral_zero( addr1(Lset(k)), addr1(Dset(m)),&
|
|
|
|
addr2(Lset(k)), addr2(Dset(m)) ) ) then
|
2024-06-04 19:31:39 +02:00
|
|
|
Delta_col(k) = &
|
2024-06-02 19:03:05 +02:00
|
|
|
ao_two_e_integral(addr1(Lset(k)), addr2(Lset(k)),&
|
|
|
|
addr1(Dset(m)), addr2(Dset(m)))
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
2024-06-04 19:31:39 +02:00
|
|
|
else
|
2024-06-17 14:57:48 +02:00
|
|
|
PROVIDE ao_integrals_map
|
2024-06-04 19:31:39 +02:00
|
|
|
!$OMP PARALLEL DO PRIVATE(k) SCHEDULE(dynamic,21)
|
|
|
|
do k=1,np
|
2024-06-17 14:57:48 +02:00
|
|
|
Delta_col(k) = 0.d0
|
2024-06-02 19:03:05 +02:00
|
|
|
if (.not.ao_two_e_integral_zero( addr1(Lset(k)), addr1(Dset(m)),&
|
|
|
|
addr2(Lset(k)), addr2(Dset(m)) ) ) then
|
2024-06-04 19:31:39 +02:00
|
|
|
Delta_col(k) = &
|
2024-06-02 19:03:05 +02:00
|
|
|
get_ao_two_e_integral( addr1(Lset(k)), addr1(Dset(m)),&
|
|
|
|
addr2(Lset(k)), addr2(Dset(m)), ao_integrals_map)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
endif
|
2024-06-04 19:31:39 +02:00
|
|
|
|
|
|
|
!$OMP PARALLEL DO PRIVATE(p)
|
|
|
|
do p=1,np
|
|
|
|
Ltmp_p(p,iblock) = Ltmp_p(p,iblock) + Delta_col(p)
|
|
|
|
Delta(p,dj) = Ltmp_p(p,iblock)
|
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
|
2023-07-11 22:31:51 +02:00
|
|
|
computed(dj) = .True.
|
|
|
|
endif
|
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! iv.
|
|
|
|
if (iblock > 1) then
|
|
|
|
call dgemv('N', np, iblock-1, -1.d0, Ltmp_p, np, Ltmp_q(dj,1), nq, 1.d0,&
|
|
|
|
Ltmp_p(1,iblock), 1)
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! iii.
|
|
|
|
f = 1.d0/dsqrt(Qmax)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-05-31 20:30:48 +02:00
|
|
|
!$OMP PARALLEL PRIVATE(p,q) DEFAULT(shared)
|
2023-07-11 17:31:58 +02:00
|
|
|
!$OMP DO
|
|
|
|
do p=1,np
|
|
|
|
Ltmp_p(p,iblock) = Ltmp_p(p,iblock) * f
|
|
|
|
L(Lset(p), rank) = Ltmp_p(p,iblock)
|
|
|
|
D(Lset(p)) = D(Lset(p)) - Ltmp_p(p,iblock) * Ltmp_p(p,iblock)
|
|
|
|
enddo
|
|
|
|
!$OMP END DO
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
!$OMP DO
|
|
|
|
do q=1,nq
|
|
|
|
Ltmp_q(q,iblock) = L(Dset(q), rank)
|
|
|
|
enddo
|
|
|
|
!$OMP END DO
|
|
|
|
!$OMP END PARALLEL
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
Qmax = D(Dset(1))
|
|
|
|
do q=1,nq
|
|
|
|
Qmax = max(Qmax, D(Dset(q)))
|
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
print '(I10, 4X, ES12.3)', rank, Qmax
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 22:31:51 +02:00
|
|
|
deallocate(Ltmp_p)
|
|
|
|
deallocate(Ltmp_q)
|
2024-06-17 14:57:48 +02:00
|
|
|
deallocate(Delta)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! i.
|
|
|
|
N = rank
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
! j.
|
2024-06-17 14:57:48 +02:00
|
|
|
D_sorted(:) = -D(:)
|
|
|
|
call dsort_noidx_big(D_sorted,ndim8)
|
|
|
|
D_sorted(:) = -D_sorted(:)
|
|
|
|
|
|
|
|
Dmax = D_sorted(1)
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-04 19:31:39 +02:00
|
|
|
np8=0_8
|
|
|
|
do p8=1,ndim8
|
2024-06-20 13:43:46 +02:00
|
|
|
if ( Dmax*D(p8) >= tau2 ) then
|
2024-06-04 19:31:39 +02:00
|
|
|
np8 = np8+1_8
|
|
|
|
Lset(np8) = p8
|
|
|
|
endif
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2024-07-03 14:52:11 +02:00
|
|
|
np = int(np8,4)
|
2024-05-31 20:50:30 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-04 16:15:09 +02:00
|
|
|
|
|
|
|
print *, '============ ============='
|
|
|
|
print *, ''
|
|
|
|
|
2024-06-20 13:43:46 +02:00
|
|
|
deallocate( D, Lset, Dset, D_sorted )
|
|
|
|
deallocate( addr1, addr2, Delta_col, computed )
|
|
|
|
|
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
allocate(cholesky_ao(ao_num,ao_num,rank), stat=ierr)
|
2024-05-31 20:07:29 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (ierr /= 0) then
|
2023-07-11 22:17:40 +02:00
|
|
|
call print_memory_usage()
|
2023-07-11 17:31:58 +02:00
|
|
|
print *, irp_here, ': Allocation failed'
|
|
|
|
stop -1
|
|
|
|
endif
|
2024-06-04 16:15:09 +02:00
|
|
|
|
|
|
|
|
|
|
|
!$OMP PARALLEL DO PRIVATE(k,j)
|
2023-07-11 17:31:58 +02:00
|
|
|
do k=1,rank
|
2024-05-31 20:50:30 +02:00
|
|
|
do j=1,ao_num
|
2024-06-17 14:57:48 +02:00
|
|
|
cholesky_ao(1:ao_num,j,k) = L((j-1_8)*ao_num+1_8:1_8*j*ao_num,k)
|
2024-05-31 20:50:30 +02:00
|
|
|
enddo
|
2023-07-11 17:31:58 +02:00
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-04 19:31:39 +02:00
|
|
|
call munmap( (/ ndim8, rank_max /), 8, fd(1), c_pointer(1) )
|
2024-06-04 16:15:09 +02:00
|
|
|
|
|
|
|
cholesky_ao_num = rank
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
if (write_ao_cholesky) then
|
2024-06-07 16:09:53 +02:00
|
|
|
print *, 'Writing Cholesky AO vectors to disk...'
|
2023-07-11 17:31:58 +02:00
|
|
|
iunit = getUnitAndOpen(trim(ezfio_work_dir)//'cholesky_ao', 'W')
|
|
|
|
write(iunit) rank
|
|
|
|
write(iunit) cholesky_ao
|
|
|
|
close(iunit)
|
|
|
|
call ezfio_set_ao_two_e_ints_io_ao_cholesky('Read')
|
|
|
|
endif
|
|
|
|
|
|
|
|
endif
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2023-07-11 17:31:58 +02:00
|
|
|
print *, 'Rank : ', cholesky_ao_num, '(', 100.d0*dble(cholesky_ao_num)/dble(ao_num*ao_num), ' %)'
|
|
|
|
print *, ''
|
2024-06-02 19:16:56 +02:00
|
|
|
call wall_time(wall1)
|
2024-06-07 16:09:53 +02:00
|
|
|
print*,'Time to provide AO cholesky vectors = ',(wall1-wall0)/60.d0, ' min'
|
2023-08-21 09:56:07 +02:00
|
|
|
|
2024-06-17 14:57:48 +02:00
|
|
|
|
2023-07-04 22:17:31 +02:00
|
|
|
END_PROVIDER
|
2023-08-21 09:56:07 +02:00
|
|
|
|