10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-11-08 15:13:48 +01:00
QuantumPackage/src/ao_two_e_ints/cholesky.irp.f

328 lines
7.2 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
do k=1,ao_num
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-04 22:17:31 +02:00
BEGIN_PROVIDER [ integer, cholesky_ao_num ]
&BEGIN_PROVIDER [ double precision, cholesky_ao, (ao_num, ao_num, 1) ]
2023-07-03 17:41:34 +02:00
implicit none
2023-07-03 18:19:31 +02:00
BEGIN_DOC
2023-07-04 22:17:31 +02:00
! 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
2023-07-03 18:19:31 +02:00
END_DOC
2023-07-04 22:17:31 +02:00
integer :: rank, ndim
double precision :: tau
double precision, pointer :: L(:,:), L_old(:,:)
2023-07-03 17:41:34 +02:00
double precision, parameter :: s = 1.d-2
double precision, parameter :: dscale = 1.d0
2023-07-03 19:54:00 +02:00
double precision, allocatable :: D(:), Delta(:,:), Ltmp_p(:,:), Ltmp_q(:,:)
2023-07-04 22:17:31 +02:00
integer, allocatable :: Lset(:), Dset(:), addr(:,:), LDmap(:), DLmap(:)
integer, allocatable :: Lset_rev(:), Dset_rev(:)
2023-07-03 17:41:34 +02:00
2023-07-04 22:17:31 +02:00
integer :: i,j,k,m,p,q, qj, dj, p2, q2
integer :: N, np, nq
2023-07-03 17:41:34 +02:00
double precision :: Dmax, Dmin, Qmax, f
2023-07-03 18:19:31 +02:00
double precision, external :: get_ao_two_e_integral
2023-07-03 21:04:50 +02:00
logical, external :: ao_two_e_integral_zero
2023-07-04 22:17:31 +02:00
integer :: block_size, iblock, ierr
PROVIDE ao_two_e_integrals_in_map
deallocate(cholesky_ao)
ndim = ao_num*ao_num
tau = ao_cholesky_threshold
allocate(L(ndim,1))
print *, ''
print *, 'Cholesky decomposition of AO integrals'
print *, '======================================'
print *, ''
print *, '============ ============='
print *, ' Rank Threshold'
print *, '============ ============='
2023-07-04 10:46:05 +02:00
2023-07-04 01:22:12 +02:00
rank = 0
2023-07-03 17:41:34 +02:00
2023-07-04 01:22:12 +02:00
allocate( D(ndim), Lset(ndim), LDmap(ndim), DLmap(ndim), Dset(ndim) )
allocate( Lset_rev(ndim), Dset_rev(ndim) )
allocate( addr(3,ndim) )
2023-07-03 17:41:34 +02:00
! 1.
2023-07-03 18:19:31 +02:00
k=0
2023-07-04 01:22:12 +02:00
do j=1,ao_num
do i=1,ao_num
2023-07-03 18:19:31 +02:00
k = k+1
addr(1,k) = i
addr(2,k) = j
2023-07-04 01:22:12 +02:00
addr(3,k) = (i-1)*ao_num + j
2023-07-03 18:19:31 +02:00
enddo
enddo
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(i)
2023-07-03 17:41:34 +02:00
do i=1,ndim
2023-07-03 18:19:31 +02:00
D(i) = get_ao_two_e_integral(addr(1,i), addr(1,i), &
addr(2,i), addr(2,i), &
ao_integrals_map)
2023-07-03 17:41:34 +02:00
enddo
2023-07-03 18:19:31 +02:00
!$OMP END PARALLEL DO
2023-07-03 17:41:34 +02:00
Dmax = maxval(D)
! 2.
np=0
2023-07-04 01:22:12 +02:00
Lset_rev = 0
2023-07-03 17:41:34 +02:00
do p=1,ndim
if ( dscale*dscale*Dmax*D(p) > tau*tau ) then
np = np+1
Lset(np) = p
2023-07-04 01:22:12 +02:00
Lset_rev(p) = np
2023-07-03 17:41:34 +02:00
endif
enddo
! 3.
N = 0
2023-07-03 18:19:31 +02:00
! 4.
2023-07-03 17:41:34 +02:00
i = 0
2023-07-03 18:19:31 +02:00
! 5.
2023-07-03 17:41:34 +02:00
do while (Dmax > tau)
! a.
i = i+1
! b.
Dmin = max(s*Dmax, tau)
! c.
nq=0
2023-07-04 01:22:12 +02:00
LDmap = 0
DLmap = 0
do p=1,np
if ( D(Lset(p)) > Dmin ) then
2023-07-03 17:41:34 +02:00
nq = nq+1
2023-07-04 01:22:12 +02:00
Dset(nq) = Lset(p)
Dset_rev(Dset(nq)) = nq
LDmap(p) = nq
DLmap(nq) = p
2023-07-03 17:41:34 +02:00
endif
enddo
2023-07-03 18:19:31 +02:00
! d., e.
2023-07-04 22:17:31 +02:00
block_size = max(N,24)
L_old => L
allocate(L(ndim,rank+nq), stat=ierr)
if (ierr /= 0) then
print *, irp_here, ': allocation failed : (Delta(np,nq))'
stop -1
endif
!$OMP PARALLEL DO PRIVATE(k)
do k=1,rank
L(:,k) = L_old(:,k)
enddo
!$OMP END PARALLEL DO
deallocate(L_old)
allocate(Delta(np,nq), stat=ierr)
if (ierr /= 0) then
print *, irp_here, ': allocation failed : (Delta(np,nq))'
stop -1
endif
allocate(Ltmp_p(np,block_size), stat=ierr)
if (ierr /= 0) then
print *, irp_here, ': allocation failed : (Ltmp_p(np,block_size))'
stop -1
endif
allocate(Ltmp_q(nq,block_size), stat=ierr)
if (ierr /= 0) then
print *, irp_here, ': allocation failed : (Ltmp_q(nq,block_size))'
stop -1
endif
2023-07-04 01:22:12 +02:00
2023-07-03 21:04:50 +02:00
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(m,k,p,q)
!$OMP DO
2023-07-03 19:54:00 +02:00
do k=1,N
do p=1,np
Ltmp_p(p,k) = L(Lset(p),k)
enddo
do q=1,nq
Ltmp_q(q,k) = L(Dset(q),k)
enddo
enddo
2023-07-04 01:22:12 +02:00
!$OMP END DO
2023-07-03 21:04:50 +02:00
2023-07-04 01:22:12 +02:00
!$OMP DO SCHEDULE(dynamic,8)
2023-07-03 21:04:50 +02:00
do m=1,nq
2023-07-04 01:22:12 +02:00
2023-07-04 22:17:31 +02:00
Delta(:,m) = 0.d0
2023-07-04 01:22:12 +02:00
do k=1, nq
! Apply only to (k,m) pairs both in Dset
p = DLmap(k)
q = Lset_rev(addr(3,Dset(k)))
if ((0 < q).and.(q < p)) cycle
2023-07-04 22:17:31 +02:00
if (.not.ao_two_e_integral_zero( addr(1,Dset(k)), addr(1,Dset(m)), &
2023-07-04 01:22:12 +02:00
addr(2,Dset(k)), addr(2,Dset(m)) ) ) then
Delta(p,m) = get_ao_two_e_integral( addr(1,Dset(k)), addr(1,Dset(m)), &
addr(2,Dset(k)), addr(2,Dset(m)), ao_integrals_map)
2023-07-04 22:17:31 +02:00
Delta(q,m) = Delta(p,m)
2023-07-04 01:22:12 +02:00
endif
2023-07-03 21:04:50 +02:00
enddo
2023-07-03 19:54:00 +02:00
2023-07-03 17:41:34 +02:00
do k=1,np
2023-07-04 01:22:12 +02:00
! Apply only to (k,m) pairs where k is not in Dset
if (LDmap(k) /= 0) cycle
q = Lset_rev(addr(3,Lset(k)))
if ((0 < q).and.(q < k)) cycle
2023-07-04 22:17:31 +02:00
if (.not.ao_two_e_integral_zero( addr(1,Lset(k)), addr(1,Dset(m)), &
2023-07-04 01:22:12 +02:00
addr(2,Lset(k)), addr(2,Dset(m)) ) ) then
Delta(k,m) = get_ao_two_e_integral( addr(1,Lset(k)), addr(1,Dset(m)), &
addr(2,Lset(k)), addr(2,Dset(m)), ao_integrals_map)
2023-07-04 22:17:31 +02:00
Delta(q,m) = Delta(k,m)
2023-07-04 01:22:12 +02:00
endif
2023-07-03 17:41:34 +02:00
enddo
enddo
2023-07-03 21:04:50 +02:00
!$OMP END DO
!$OMP END PARALLEL
2023-07-03 17:41:34 +02:00
2023-07-04 22:17:31 +02:00
if (N>0) then
call dgemm('N','T', np, nq, N, -1.d0, &
Ltmp_p, np, Ltmp_q, nq, 1.d0, Delta, np)
endif
2023-07-03 19:54:00 +02:00
! f.
2023-07-03 17:41:34 +02:00
Qmax = D(Dset(1))
do q=1,nq
Qmax = max(Qmax, D(Dset(q)))
enddo
! g.
2023-07-04 10:46:05 +02:00
iblock = 0
2023-07-04 01:22:12 +02:00
do j=1,nq
if ( (Qmax <= Dmin).or.(N+j > ndim) ) exit
2023-07-03 17:41:34 +02:00
! i.
rank = N+j
2023-07-04 10:46:05 +02:00
if (iblock == block_size) then
call dgemm('N','T',np,nq,block_size,-1.d0, &
Ltmp_p, np, Ltmp_q, nq, 1.d0, Delta, np)
iblock = 0
endif
2023-07-03 17:41:34 +02:00
! ii.
do dj=1,nq
qj = Dset(dj)
if (D(qj) == Qmax) then
exit
endif
enddo
2023-07-04 22:17:31 +02:00
L(1:ndim, rank) = 0.d0
2023-07-04 01:22:12 +02:00
2023-07-04 10:46:05 +02:00
iblock = iblock+1
do p=1,np
Ltmp_p(p,iblock) = Delta(p,dj)
enddo
2023-07-04 22:17:31 +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, &
2023-07-04 10:46:05 +02:00
Ltmp_p(1,iblock), 1)
2023-07-04 22:17:31 +02:00
endif
2023-07-04 10:46:05 +02:00
2023-07-03 17:41:34 +02:00
! iii.
f = 1.d0/dsqrt(Qmax)
2023-07-04 01:22:12 +02:00
!$OMP PARALLEL PRIVATE(m,p,q,k) DEFAULT(shared)
2023-07-03 21:04:50 +02:00
!$OMP DO
2023-07-03 17:41:34 +02:00
do p=1,np
2023-07-04 10:46:05 +02:00
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)
2023-07-03 19:54:00 +02:00
enddo
2023-07-03 21:04:50 +02:00
!$OMP END DO
2023-07-03 19:54:00 +02:00
2023-07-03 21:04:50 +02:00
!$OMP DO
2023-07-03 19:54:00 +02:00
do q=1,nq
2023-07-04 10:46:05 +02:00
Ltmp_q(q,iblock) = L(Dset(q), rank)
2023-07-03 17:41:34 +02:00
enddo
2023-07-03 21:04:50 +02:00
!$OMP END DO
2023-07-04 01:22:12 +02:00
2023-07-03 21:04:50 +02:00
!$OMP END PARALLEL
2023-07-03 17:41:34 +02:00
Qmax = D(Dset(1))
2023-07-04 22:17:31 +02:00
do q=1,nq
Qmax = max(Qmax, D(Dset(q)))
2023-07-03 17:41:34 +02:00
enddo
enddo
2023-07-04 22:17:31 +02:00
print '(I10, 4X, ES12.3)', rank, Qmax
deallocate(Delta, stat=ierr)
deallocate(Ltmp_p, stat=ierr)
deallocate(Ltmp_q, stat=ierr)
2023-07-03 17:41:34 +02:00
! i.
N = N+j
! j.
Dmax = D(Lset(1))
do p=1,np
Dmax = max(Dmax, D(Lset(p)))
enddo
np=0
2023-07-04 01:22:12 +02:00
Lset_rev = 0
2023-07-03 17:41:34 +02:00
do p=1,ndim
if ( dscale*dscale*Dmax*D(p) > tau*tau ) then
np = np+1
Lset(np) = p
2023-07-04 01:22:12 +02:00
Lset_rev(p) = np
2023-07-03 17:41:34 +02:00
endif
enddo
2023-07-03 18:19:31 +02:00
2023-07-03 17:41:34 +02:00
enddo
2023-07-04 22:17:31 +02:00
allocate(cholesky_ao(ao_num,ao_num,rank))
call dcopy(ndim*rank, L, 1, cholesky_ao, 1)
deallocate(L)
cholesky_ao_num = rank
print *, '============ ============='
print *, ''
print *, 'Rank : ', cholesky_ao_num, '(', 100.d0*dble(cholesky_ao_num)/dble(ao_num*ao_num), ' %)'
print *, ''
END_PROVIDER