4
1
mirror of https://github.com/pfloos/quack synced 2024-11-12 17:13:57 +01:00
quack/src/utils/orthogonalization_matrix.f90

119 lines
2.1 KiB
Fortran
Raw Normal View History

2019-03-19 10:13:33 +01:00
subroutine orthogonalization_matrix(ortho_type,nBas,S,X)
! Compute the orthogonalization matrix X
implicit none
! Input variables
integer,intent(in) :: nBas,ortho_type
double precision,intent(in) :: S(nBas,nBas)
! Local variables
logical :: debug
double precision,allocatable :: UVec(:,:),Uval(:)
2020-06-01 11:35:17 +02:00
double precision,parameter :: thresh = 1d-6
2019-03-19 10:13:33 +01:00
integer :: i
! Output variables
double precision,intent(out) :: X(nBas,nBas)
debug = .false.
! Type of orthogonalization ortho_type
!
! 1 = Lowdin
! 2 = Canonical
! 3 = SVD
!
allocate(Uvec(nBas,nBas),Uval(nBas))
if(ortho_type == 1) then
2020-03-25 21:10:21 +01:00
! write(*,*)
! write(*,*) ' Lowdin orthogonalization'
! write(*,*)
2019-03-19 10:13:33 +01:00
Uvec = S
call diagonalize_matrix(nBas,Uvec,Uval)
do i=1,nBas
2020-05-21 16:13:52 +02:00
if(Uval(i) < thresh) then
2019-03-19 10:13:33 +01:00
2020-06-01 11:35:17 +02:00
write(*,*) 'Eigenvalue',i,' is very small in Lowdin orthogonalization = ',Uval(i)
2019-03-19 10:13:33 +01:00
endif
2020-05-21 16:13:52 +02:00
Uval(i) = 1d0/sqrt(Uval(i))
2019-03-19 10:13:33 +01:00
enddo
call ADAt(nBas,Uvec,Uval,X)
elseif(ortho_type == 2) then
2021-02-16 10:00:01 +01:00
! write(*,*)
! write(*,*) 'Canonical orthogonalization'
! write(*,*)
2019-03-19 10:13:33 +01:00
Uvec = S
call diagonalize_matrix(nBas,Uvec,Uval)
do i=1,nBas
if(Uval(i) > thresh) then
Uval(i) = 1d0/sqrt(Uval(i))
else
write(*,*) ' Eigenvalue',i,'too small for canonical orthogonalization'
endif
enddo
call AD(nBas,Uvec,Uval)
X = Uvec
elseif(ortho_type == 3) then
2021-02-16 10:00:01 +01:00
! write(*,*)
! write(*,*) ' SVD-based orthogonalization NYI'
! write(*,*)
2019-03-19 10:13:33 +01:00
! Uvec = S
! call diagonalize_matrix(nBas,Uvec,Uval)
! do i=1,nBas
! if(Uval(i) > thresh) then
! Uval(i) = 1d0/sqrt(Uval(i))
! else
! write(*,*) 'Eigenvalue',i,'too small for canonical orthogonalization'
! endif
! enddo
!
! call AD(nBas,Uvec,Uval)
! X = Uvec
endif
! Print results
if(debug) then
write(*,'(A28)') '----------------------'
write(*,'(A28)') 'Orthogonalization matrix'
write(*,'(A28)') '----------------------'
call matout(nBas,nBas,X)
write(*,*)
endif
end subroutine orthogonalization_matrix