10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-22 20:34:58 +01:00

Accelerated cholesky AO-MO transformation

This commit is contained in:
Anthony Scemama 2023-07-05 19:12:03 +02:00
parent 41a369dd68
commit 5a0c8de5a3

View File

@ -4,16 +4,18 @@ BEGIN_PROVIDER [ double precision, cholesky_mo, (mo_num, mo_num, cholesky_ao_num
! Cholesky vectors in MO basis
END_DOC
integer :: k
integer :: k, i, j
call set_multiple_levels_omp(.False.)
print *, 'AO->MO Transformation of Cholesky vectors'
!$OMP PARALLEL DO PRIVATE(k)
do k=1,cholesky_ao_num
call ao_to_mo(cholesky_ao(1,1,k),ao_num,cholesky_mo(1,1,k),mo_num)
do j=1,mo_num
do i=1,mo_num
cholesky_mo(i,j,k) = cholesky_mo_transp(k,i,j)
enddo
enddo
enddo
!$OMP END PARALLEL DO
print *, ''
END_PROVIDER
@ -23,27 +25,15 @@ BEGIN_PROVIDER [ double precision, cholesky_mo_transp, (cholesky_ao_num, mo_num,
! Cholesky vectors in MO basis
END_DOC
integer :: i,j,k
double precision, allocatable :: buffer(:,:)
double precision, allocatable :: X(:,:,:)
print *, 'AO->MO Transformation of Cholesky vectors'
print *, 'AO->MO Transformation of Cholesky vectors .'
call set_multiple_levels_omp(.False.)
!$OMP PARALLEL PRIVATE(i,j,k,buffer)
allocate(buffer(mo_num,mo_num))
!$OMP DO SCHEDULE(static)
do k=1,cholesky_ao_num
call ao_to_mo(cholesky_ao(1,1,k),ao_num,buffer,mo_num)
do j=1,mo_num
do i=1,mo_num
cholesky_mo_transp(k,i,j) = buffer(i,j)
enddo
enddo
enddo
!$OMP END DO
deallocate(buffer)
!$OMP END PARALLEL
print *, ''
allocate(X(mo_num,cholesky_ao_num,ao_num))
call dgemm('T','N', ao_num*cholesky_ao_num, mo_num, ao_num, 1.d0, &
cholesky_ao, ao_num, mo_coef, ao_num, 0.d0, X, ao_num*cholesky_ao_num)
call dgemm('T','N', cholesky_ao_num*mo_num, mo_num, ao_num, 1.d0, &
X, ao_num, mo_coef, ao_num, 0.d0, cholesky_mo_transp, cholesky_ao_num*mo_num)
deallocate(X)
END_PROVIDER