10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-11-15 10:33:50 +01:00

Rewrote Cholesky

This commit is contained in:
Anthony Scemama 2023-07-03 17:41:34 +02:00
parent f30ed6bf43
commit d911f4eee8

View File

@ -1,10 +1,28 @@
BEGIN_PROVIDER [ integer, cholesky_ao_num_guess ]
BEGIN_PROVIDER [ integer, mini_basis_size, (128) ]
implicit none
BEGIN_DOC
! Size of the minimal basis set per element
END_DOC
mini_basis_size(1:2) = 1
mini_basis_size(3:4) = 2
mini_basis_size(5:10) = 5
mini_basis_size(11:12) = 6
mini_basis_size(13:18) = 9
mini_basis_size(19:20) = 13
mini_basis_size(21:36) = 18
mini_basis_size(37:38) = 22
mini_basis_size(39:54) = 27
mini_basis_size(55:) = 36
END_PROVIDER
BEGIN_PROVIDER [ integer, cholesky_ao_num_guess ]
implicit none
BEGIN_DOC
! Number of Cholesky vectors in AO basis
END_DOC
cholesky_ao_num_guess = ao_num*ao_num
cholesky_ao_num_guess = ao_num*ao_num !sum(mini_basis_size(int(nucl_charge(:))))
END_PROVIDER
BEGIN_PROVIDER [ integer, cholesky_ao_num ]
@ -103,8 +121,10 @@ END_PROVIDER
! Call Lapack
cholesky_ao_num = cholesky_ao_num_guess
call pivoted_cholesky(ao_integrals, cholesky_ao_num, ao_cholesky_threshold, ao_num*ao_num, cholesky_ao)
print *, 'Rank: ', cholesky_ao_num, '(', 100.d0*dble(cholesky_ao_num)/dble(ao_num*ao_num), ' %)'
! call pivoted_cholesky(ao_integrals, cholesky_ao_num, ao_cholesky_threshold, ao_num*ao_num, cholesky_ao)
call direct_cholesky(ao_integrals, cholesky_ao_num, ao_cholesky_threshold, ao_num*ao_num, cholesky_ao)
print *, 'Rank : ', cholesky_ao_num, '(', 100.d0*dble(cholesky_ao_num)/dble(ao_num*ao_num), ' %)'
! Remove mmap
double precision, external :: getUnitAndOpen
@ -131,3 +151,172 @@ BEGIN_PROVIDER [ double precision, cholesky_ao_transp, (cholesky_ao_num, ao_num,
enddo
END_PROVIDER
subroutine direct_cholesky( A, rank, tau, ndim, L)
implicit none
integer :: ndim
integer, intent(inout) :: rank
double precision, intent(inout) :: A(ndim, ndim)
double precision, intent(out) :: L(ndim, rank)
double precision, intent(in) :: tau
double precision, parameter :: s = 1.d-2
double precision, parameter :: dscale = 1.d0
double precision, allocatable :: D(:), Delta(:,:)
integer, allocatable :: Lset(:), Dset(:)
integer :: i,j,k,m,p,q, qj, dj
integer :: N, np, nq
double precision :: Dmax, Dmin, Qmax, f
allocate( D(ndim), Lset(ndim), Dset(ndim) )
L = 0.d0
! 1.
do i=1,ndim
D(i) = A(i,i)
enddo
Dmax = maxval(D)
! print *, '# 1. ', D
! print *, '# 1. ', Dmax
! 2.
np=0
do p=1,ndim
if ( dscale*dscale*Dmax*D(p) > tau*tau ) then
np = np+1
Lset(np) = p
endif
enddo
! print *, '# 2. ', Lset(:np)
! 3.
N = 0
! print *, '# 3. ', N
! 4.
i = 0
! print *, '# 4. ', i
! 5.
do while (Dmax > tau)
! a.
i = i+1
! print *, '# 5.a ', i
! b.
Dmin = max(s*Dmax, tau)
! print *, '# 5.b ', Dmin
! c.
nq=0
do q=1,np
if ( D(Lset(q)) > Dmin ) then
nq = nq+1
Dset(nq) = Lset(q)
endif
enddo
! print *, '# 5.c ', Dset(:nq)
! d.
allocate(Delta(np,nq))
do m=1,nq
do k=1,np
Delta(k,m) = A(Lset(k), Dset(m))
enddo
enddo
! print *, '# 5.d ', Delta
! e.
do m=1,nq
do k=1,np
do p=1,N
Delta(k,m) = Delta(k,m) - L(Lset(k),p) * L(Dset(m),p)
enddo
enddo
enddo
! print *, '# 5.e ', Delta
! f.
Qmax = D(Dset(1))
do q=1,nq
Qmax = max(Qmax, D(Dset(q)))
enddo
! print *, '# 5.f ', Qmax
! g.
j = 0
! print *, '# 5.g ', j
do while ( (j <= nq).and.(Qmax > Dmin) )
! i.
j = j+1
rank = N+j
! print *, '# 5.h.i ', j, rank
! ii.
do dj=1,nq
qj = Dset(dj)
if (D(qj) == Qmax) then
exit
endif
enddo
! print *, ' # 5.h.ii ', qj, dj
! iii.
f = 1.d0/dsqrt(Qmax)
do p=1,np
L(Lset(p), rank) = Delta(p,dj) * f
enddo
! print *, ' # 5.h.iii '
! do k=1,20
! print *, L(k,1:rank)
! enddo
! iv.
do m=1, nq
do k=1, np
Delta(k,m) = Delta(k,m) - L(Lset(k),rank) * L(Dset(m),rank)
enddo
enddo
do k=1, np
D(Lset(k)) = D(Lset(k)) - L(Lset(k),rank) * L(Lset(k),rank)
enddo
Qmax = D(Dset(1))
do q=1,np
Qmax = max(Qmax, D(Lset(q)))
enddo
! print *, '# 5.h.iv ', Delta
! print *, '# 5.h.iv ', D
! print *, '# 5.h.iv ', Qmax
enddo
deallocate(Delta)
! i.
N = N+j
! print *, '# 5.i ', N
! j.
Dmax = D(Lset(1))
do p=1,np
Dmax = max(Dmax, D(Lset(p)))
enddo
! print *, '# 5.j ', Dmax
np=0
do p=1,ndim
if ( dscale*dscale*Dmax*D(p) > tau*tau ) then
np = np+1
Lset(np) = p
endif
enddo
! print *, '# k. ', Lset(:np)
enddo
end