10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-11-20 04:52:36 +01:00
QuantumPackage/src/ao_two_e_ints/cholesky.irp.f

499 lines
14 KiB
Fortran
Raw Normal View History

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
!
! https://mogp-emulator.readthedocs.io/en/latest/methods/proc/ProcPivotedCholesky.html
! https://doi.org/10.1016/j.apnum.2011.10.001 : Page 4, Algorithm 1
2023-07-11 17:31:58 +02:00
END_DOC
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-04 19:31:39 +02:00
double precision, pointer :: L(:,:), Delta(:,:)
2023-07-11 17:31:58 +02:00
double precision :: s
2024-06-02 19:03:05 +02:00
double precision :: dscale, dscale_tmp
2024-06-04 19:31:39 +02:00
double precision, allocatable :: D(:), Ltmp_p(:,:), Ltmp_q(:,:), D_sorted(:), Delta_col(:)
2024-06-02 19:03:05 +02:00
integer, allocatable :: addr1(:), addr2(:)
integer*8, allocatable :: Lset(:), Dset(:), addr3(:)
2023-07-11 22:31:51 +02:00
logical, allocatable :: computed(:)
2024-05-31 20:50:30 +02:00
integer :: i,j,k,m,p,q, dj, p2, q2
2024-06-04 16:15:09 +02:00
integer*8 :: i8, j8, p8, qj8, rank_max, np8
integer :: N, np, nq
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-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
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-07-11 17:31:58 +02:00
integer, external :: getUnitAndOpen
2024-05-31 20:50:30 +02:00
integer :: iunit, ierr
2024-05-31 20:50:30 +02:00
ndim8 = ao_num*ao_num*1_8
double precision :: wall0,wall1
2024-06-04 16:15:09 +02:00
type(c_ptr) :: c_pointer(2)
integer :: fd(2)
2024-06-04 19:31:39 +02:00
logical :: delta_on_disk
2024-06-04 16:15:09 +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-07-11 17:31:58 +02:00
if (read_ao_cholesky) then
print *, 'Reading Cholesky vectors from disk...'
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-07-11 17:31:58 +02:00
else
2024-03-20 09:20:21 +01:00
PROVIDE nucl_coord ao_two_e_integral_schwartz
call set_multiple_levels_omp(.False.)
2024-06-04 19:31:39 +02:00
call resident_memory(mem0)
2024-06-04 16:15:09 +02:00
rank_max = min(ndim8,274877906944_8/1_8/ndim8)
call mmap(trim(ezfio_work_dir)//'cholesky_ao_tmp', (/ ndim8, rank_max /), 8, fd(1), .False., 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-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-07-11 17:31:58 +02:00
tau = ao_cholesky_threshold
2024-06-02 19:03:05 +02:00
tau2 = tau*tau
2024-05-31 20:50:30 +02:00
mem = 6.d0 * memory_of_double8(ndim8) + 6.d0 * memory_of_int8(ndim8)
2023-07-11 22:17:40 +02:00
call check_mem(mem, irp_here)
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-07-11 17:31:58 +02:00
rank = 0
2024-06-04 16:15:09 +02:00
allocate( D(ndim8), Lset(ndim8), Dset(ndim8), D_sorted(ndim8))
2024-06-02 19:03:05 +02:00
allocate( addr1(ndim8), addr2(ndim8), addr3(ndim8) )
!print *, 'allocate : (D(ndim8))', memory_of_int8(ndim8)
!print *, 'allocate : (Lset(ndim8))', memory_of_int8(ndim8)
!print *, 'allocate : (Dset(ndim8))', memory_of_int8(ndim8)
!print *, 'allocate : (4,addr(ndim8))', memory_of_int8(4_8*ndim8)
2023-07-11 17:31:58 +02:00
! 1.
k=0
do j=1,ao_num
do i=1,ao_num
k = k+1
2024-05-31 20:50:30 +02:00
addr1(k) = i
addr2(k) = j
addr3(k) = (i-1)*ao_num + j
2023-07-11 17:31:58 +02:00
enddo
enddo
2023-07-11 17:31:58 +02:00
if (do_direct_integrals) then
2024-06-02 19:03:05 +02:00
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(i8) SCHEDULE(dynamic,16)
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-02 19:03:05 +02:00
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(i8) SCHEDULE(dynamic,16)
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), &
addr2(i8), addr2(i8), &
2023-07-11 17:31:58 +02:00
ao_integrals_map)
enddo
!$OMP END PARALLEL DO
endif
2024-06-04 16:15:09 +02:00
D_sorted(:) = -D(:)
call dsort_noidx_big(D_sorted,ndim8)
D_sorted(:) = dabs(D_sorted(:))
2024-06-04 16:15:09 +02:00
Dmax = D_sorted(1)
2023-07-11 17:31:58 +02:00
! 2.
2024-06-04 19:31:39 +02:00
dscale = 1.d0
dscale_tmp = dscale*dscale*Dmax
np8=0_8
do p8=1,ndim8
if ( dscale_tmp*D(p8) > tau2 ) then
np8 = np8+1_8
Lset(np8) = p8
endif
2023-07-11 17:31:58 +02:00
enddo
2024-06-04 19:31:39 +02:00
np = np8
2023-07-11 17:31:58 +02:00
! 3.
N = 0
2023-07-11 17:31:58 +02:00
! 4.
i = 0
2023-07-11 17:31:58 +02:00
! 5.
2024-06-04 16:15:09 +02:00
do while ( (Dmax > tau).and.(rank*1_8 < min(ndim8,rank_max)) )
2023-07-11 17:31:58 +02:00
! a.
i = i+1
2023-07-11 17:31:58 +02:00
! Inrease s until the arrays fit in memory
s = 0.01d0
block_size = max(N,24)
2023-07-11 17:31:58 +02:00
do while (.True.)
2023-07-11 17:31:58 +02:00
! b.
Dmin = max(s*Dmax,tau)
2023-07-11 17:31:58 +02:00
! c.
nq=0
do p=1,np
if ( D(Lset(p)) > Dmin ) then
nq = nq+1
Dset(nq) = Lset(p)
2023-07-06 17:51:59 +02:00
endif
2023-07-11 17:31:58 +02:00
enddo
2024-06-04 19:31:39 +02:00
mem = mem0 &
+ np*memory_of_double(nq)
2024-06-04 16:15:09 +02:00
!print *, 'mem = ', mem
2024-06-04 19:31:39 +02:00
if (mem > 300.d0) then ! 300GB max for Delta
2023-07-11 17:31:58 +02:00
s = s*2.d0
else
exit
endif
2023-07-11 17:31:58 +02:00
if ((s > 1.d0).or.(nq == 0)) then
2023-07-11 22:17:40 +02:00
call print_memory_usage()
print *, 'Required peak memory: ', mem, 'Gb'
2024-06-04 16:15:09 +02:00
call resident_memory(mem)
print *, 'Already used memory: ', mem, 'Gb'
print *, 'Not enough memory. Reduce cholesky threshold'
2023-07-11 17:31:58 +02:00
stop -1
endif
2024-06-04 19:31:39 +02:00
if (s > 0.1d0) then
exit
endif
2023-07-11 17:31:58 +02:00
enddo
2023-07-11 17:31:58 +02:00
! d., e.
2024-06-04 19:31:39 +02:00
mem = mem0 &
+ memory_of_int(nq) &! computed(nq)
+ np*memory_of_int(nq) &! computed(nq)
+ (np+nq)*memory_of_double(block_size) ! Ltmp_p(np,block_size) + Ltmp_q(nq,block_size)
2024-06-04 19:31:39 +02:00
if (mem > qp_max_mem) then
call mmap(trim(ezfio_work_dir)//'cholesky_delta', (/ np*1_8, nq*1_8 /), 8, fd(2), .False., c_pointer(2))
call c_f_pointer(c_pointer(2), Delta, (/ np, nq /))
! 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_delta', 'R')
close(iunit,status='delete')
delta_on_disk = .True.
else
allocate(Delta(np,nq))
delta_on_disk = .False.
2023-07-11 17:31:58 +02:00
endif
2024-06-05 02:51:12 +02:00
!print *, delta_on_disk
2024-06-04 19:31:39 +02:00
allocate(Delta_col(np))
2023-07-11 17:31:58 +02:00
allocate(Ltmp_p(np,block_size), stat=ierr)
2024-06-02 19:03:05 +02:00
!print *, 'allocate : Ltmp_p(np,block_size)', memory_of_double8(np*block_size*1_8), np, block_size
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-07-11 17:31:58 +02:00
allocate(Ltmp_q(nq,block_size), stat=ierr)
2024-06-02 19:03:05 +02:00
!print *, 'allocate : Ltmp_q(nq,block_size)', memory_of_double8(nq*block_size*1_8), nq, block_size
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-07-11 22:31:51 +02:00
allocate(computed(nq))
2024-06-05 03:16:55 +02:00
computed(:) = .False.
2024-06-05 03:16:55 +02:00
!print *, 'N, rank, block_size', N, rank, block_size
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-07-11 17:31:58 +02:00
if (N>0) then
2024-06-04 19:31:39 +02:00
call dgemm('N','T', np, nq, N, -1.d0, &
Ltmp_p, np, Ltmp_q, nq, 0.d0, Delta, np)
else
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(q,j) SCHEDULE(static,1)
do q=1,nq
do j=1,np
Delta(j,q) = 0.d0
enddo
enddo
!$OMP END PARALLEL DO
2023-07-11 17:31:58 +02:00
endif
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-07-11 17:31:58 +02:00
! g.
2023-07-11 17:31:58 +02:00
iblock = 0
do j=1,nq
2024-05-31 20:50:30 +02:00
if ( (Qmax <= Dmin).or.(N+j*1_8 > ndim8) ) exit
2023-07-11 17:31:58 +02:00
! i.
rank = N+j
2023-07-11 17:31:58 +02:00
if (iblock == block_size) then
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)
2023-07-11 17:31:58 +02:00
iblock = 0
endif
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
2024-06-02 19:03:05 +02:00
do i8=1,ndim8
L(i8, rank) = 0.d0
enddo
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-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)))
2024-06-05 02:51:12 +02:00
else
Delta_col(k) = 0.d0
2024-06-02 19:03:05 +02:00
endif
enddo
!$OMP END PARALLEL DO
2024-06-04 19:31:39 +02:00
else
!$OMP PARALLEL DO PRIVATE(k) SCHEDULE(dynamic,21)
do k=1,np
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)
2024-06-05 02:51:12 +02:00
else
Delta_col(k) = 0.d0
2024-06-02 19:03:05 +02:00
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-07-11 17:31:58 +02:00
! iii.
f = 1.d0/dsqrt(Qmax)
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-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-07-11 17:31:58 +02:00
Qmax = D(Dset(1))
do q=1,nq
Qmax = max(Qmax, D(Dset(q)))
enddo
2023-07-11 17:31:58 +02:00
enddo
2023-07-11 17:31:58 +02:00
print '(I10, 4X, ES12.3)', rank, Qmax
2024-06-04 19:31:39 +02:00
deallocate(Delta_col)
2023-07-11 22:31:51 +02:00
deallocate(Ltmp_p)
deallocate(Ltmp_q)
2024-06-02 19:03:05 +02:00
deallocate(computed)
2024-06-04 19:31:39 +02:00
if (delta_on_disk) then
call munmap( (/ np*1_8, nq*1_8 /), 8, fd(2), c_pointer(2) )
else
deallocate(Delta)
endif
2023-07-11 17:31:58 +02:00
! i.
N = rank
2023-07-11 17:31:58 +02:00
! j.
Dmax = D(Lset(1))
do p=1,np
Dmax = max(Dmax, D(Lset(p)))
enddo
2024-06-04 19:31:39 +02:00
dscale = 1.d0
dscale_tmp = dscale*dscale*Dmax
np8=0_8
do p8=1,ndim8
if ( dscale_tmp*D(p8) > tau2 ) then
np8 = np8+1_8
Lset(np8) = p8
endif
2023-07-11 17:31:58 +02:00
enddo
2024-06-04 19:31:39 +02:00
np = np8
2024-05-31 20:50:30 +02:00
2023-07-11 17:31:58 +02:00
enddo
2024-06-04 16:15:09 +02:00
print *, '============ ============='
print *, ''
2023-07-11 17:31:58 +02:00
allocate(cholesky_ao(ao_num,ao_num,rank), stat=ierr)
2024-06-02 19:03:05 +02:00
!print *, 'allocate : cholesky_ao(ao_num,ao_num,rank)', memory_of_double8(ao_num*ao_num*rank*1_8)
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-02 19:03:05 +02:00
cholesky_ao(1:ao_num,j,k) = L((j-1)*ao_num+1: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
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-07-11 17:31:58 +02:00
if (write_ao_cholesky) then
print *, 'Writing Cholesky vectors to disk...'
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-07-11 17:31:58 +02:00
print *, 'Rank : ', cholesky_ao_num, '(', 100.d0*dble(cholesky_ao_num)/dble(ao_num*ao_num), ' %)'
print *, ''
call wall_time(wall1)
print*,'Time to provide AO cholesky vectors = ',wall1-wall0
2023-07-04 22:17:31 +02:00
END_PROVIDER