4
1
mirror of https://github.com/pfloos/quack synced 2024-11-07 06:33:55 +01:00
quack/src/utils/read_basis.f90

186 lines
5.2 KiB
Fortran
Raw Normal View History

2020-03-25 09:48:58 +01:00
subroutine read_basis(nNuc,rNuc,nBas,nO,nV,nShell,TotAngMomShell,CenterShell,KShell,DShell,ExpShell, &
max_ang_mom,min_exponent,max_exponent)
2019-03-20 13:38:42 +01:00
! Read basis set information
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: nNuc,nO
double precision,intent(in) :: rNuc(nNuc,ncart)
! Local variables
integer :: nShAt,iNuc,iShell
2020-01-17 13:45:02 +01:00
integer :: i,j,k,kk
2019-03-20 13:38:42 +01:00
character :: shelltype
! Output variables
integer,intent(out) :: nShell,nBas
double precision,intent(out) :: CenterShell(maxShell,ncart)
integer,intent(out) :: TotAngMomShell(maxShell),KShell(maxShell)
double precision,intent(out) :: DShell(maxShell,maxK),ExpShell(maxShell,maxK)
integer,intent(out) :: nV
2020-03-25 09:48:58 +01:00
integer,intent(out) :: max_ang_mom(nNuc)
double precision,intent(out) :: min_exponent(nNuc,maxL+1)
double precision,intent(out) :: max_exponent(nNuc)
2019-03-20 13:38:42 +01:00
!------------------------------------------------------------------------
! Primary basis set information
!------------------------------------------------------------------------
! Open file with basis set specification
open(unit=2,file='input/basis')
! Read basis information
write(*,'(A28)') 'Gaussian basis set'
write(*,'(A28)') '------------------'
2020-03-25 09:48:58 +01:00
! Initailization
2019-03-20 13:38:42 +01:00
nShell = 0
2020-03-25 09:48:58 +01:00
max_ang_mom(:) = 0
min_exponent(:,:) = huge(0d0)
max_exponent(:) = 0d0
!------------------------------------------------------------------------
! Loop over atoms
!------------------------------------------------------------------------
2019-03-20 13:38:42 +01:00
do i=1,nNuc
2020-03-25 09:48:58 +01:00
2019-03-20 13:38:42 +01:00
read(2,*) iNuc,nShAt
write(*,'(A28,1X,I16)') 'Atom n. ',iNuc
write(*,'(A28,1X,I16)') 'number of shells ',nShAt
write(*,'(A28)') '------------------'
2020-03-25 09:48:58 +01:00
!------------------------------------------------------------------------
! Loop over shells
!------------------------------------------------------------------------
2019-03-20 13:38:42 +01:00
do j=1,nShAt
2020-03-25 09:48:58 +01:00
2019-03-20 13:38:42 +01:00
nShell = nShell + 1
2020-03-25 09:48:58 +01:00
! Basis function centers
2019-03-20 13:38:42 +01:00
do k=1,ncart
CenterShell(nShell,k) = rNuc(iNuc,k)
enddo
2020-03-25 09:48:58 +01:00
! Shell type and contraction degree
2019-03-20 13:38:42 +01:00
read(2,*) shelltype,KShell(nShell)
2020-03-25 09:48:58 +01:00
select case (shelltype)
case ("S")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 0
write(*,'(A28,1X,I16)') 's-type shell with K = ',KShell(nShell)
2020-03-25 09:48:58 +01:00
case ("P")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 1
2020-03-25 09:48:58 +01:00
write(*,'(A28,1X,I16)') 'p-type shell with K = ',KShell(nShell)
case ("D")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 2
write(*,'(A28,1X,I16)') 'd-type shell with K = ',KShell(nShell)
2020-03-25 09:48:58 +01:00
case ("F")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 3
write(*,'(A28,1X,I16)') 'f-type shell with K = ',KShell(nShell)
2020-03-25 09:48:58 +01:00
case ("G")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 4
write(*,'(A28,1X,I16)') 'g-type shell with K = ',KShell(nShell)
2020-03-25 09:48:58 +01:00
case ("H")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 5
write(*,'(A28,1X,I16)') 'h-type shell with K = ',KShell(nShell)
2020-03-25 09:48:58 +01:00
case ("I")
2019-03-20 13:38:42 +01:00
TotAngMomShell(nShell) = 6
write(*,'(A28,1X,I16)') 'i-type shell with K = ',KShell(nShell)
2020-03-25 09:48:58 +01:00
case ("J")
TotAngMomShell(nShell) = 7
write(*,'(A28,1X,I16)') 'j-type shell with K = ',KShell(nShell)
case default
call print_warning('!!! Angular momentum too high !!!')
stop
end select
2019-03-20 13:38:42 +01:00
! Read exponents and contraction coefficients
2020-03-25 09:48:58 +01:00
write(*,'(A28,1X,A16,A16)') '','Exponents','Contraction'
do k=1,Kshell(nShell)
read(2,*) kk,ExpShell(nShell,k),DShell(nShell,k)
write(*,'(A28,1X,F16.10,F16.10)') '',ExpShell(nShell,k),DShell(nShell,k)
enddo
min_exponent(iNuc,TotAngMomShell(nShell)+1) &
= min(min_exponent(iNuc,TotAngMomShell(nShell)+1),minval(ExpShell(nShell,1:KShell(nShell))))
max_exponent(iNuc) = max(max_exponent(iNuc),maxval(ExpShell(nShell,:)))
max_ang_mom(iNuc) = max(max_ang_mom(iNuc),TotAngMomShell(nShell))
2019-03-20 13:38:42 +01:00
enddo
2020-03-25 09:48:58 +01:00
!------------------------------------------------------------------------
! End loop over shells
!------------------------------------------------------------------------
2019-03-20 13:38:42 +01:00
write(*,'(A28)') '------------------'
2020-03-25 09:48:58 +01:00
! print*,'maximum angular momemtum for atom n. ',iNuc,' = '
! print*,max_ang_mom(iNuc)
! print*,'maximum exponent for atom n. ',iNuc,' = '
! print*,max_exponent(iNuc)
! print*,'minimum exponent for atom n. ',iNuc,' = '
! print*,min_exponent(iNuc,1:max_ang_mom(iNuc)+1)
2019-03-20 13:38:42 +01:00
enddo
2020-03-25 09:48:58 +01:00
!------------------------------------------------------------------------
! End loop over atoms
!------------------------------------------------------------------------
2019-03-20 13:38:42 +01:00
! Total number of shells
write(*,'(A28,1X,I16)') 'Number of shells',nShell
write(*,'(A28)') '------------------'
write(*,*)
! Close file with basis set specification
close(unit=2)
! Calculate number of basis functions
nBas = 0
do iShell=1,nShell
nBas = nBas + (TotAngMomShell(iShell)*TotAngMomShell(iShell) + 3*TotAngMomShell(iShell) + 2)/2
enddo
write(*,'(A28)') '------------------'
write(*,'(A28,1X,I16)') 'Number of basis functions',NBas
write(*,'(A28)') '------------------'
write(*,*)
! Number of virtual orbitals
nV = nBas - nO
end subroutine read_basis