4
1
mirror of https://github.com/pfloos/quack synced 2024-11-13 09:34:04 +01:00
quack/src/LR/print_transition_vectors.f90

87 lines
2.3 KiB
Fortran
Raw Normal View History

subroutine print_transition_vectors(spin_allowed,nBas,nC,nO,nV,nR,nS,dipole_int,Omega,XpY,XmY)
2020-06-03 12:06:16 +02:00
! Print transition vectors for linear response calculation
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: spin_allowed
2020-06-03 12:06:16 +02:00
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
integer,intent(in) :: nS
double precision :: dipole_int(nBas,nBas,ncart)
2020-06-03 12:06:16 +02:00
double precision,intent(in) :: Omega(nS)
double precision,intent(in) :: XpY(nS,nS)
double precision,intent(in) :: XmY(nS,nS)
! Local variables
2020-10-05 16:58:19 +02:00
integer :: ia,jb,j,b
integer :: maxS = 10
2020-10-04 14:22:38 +02:00
double precision :: S2
double precision,parameter :: thres_vec = 0.1d0
2020-06-03 12:06:16 +02:00
double precision,allocatable :: X(:)
double precision,allocatable :: Y(:)
double precision,allocatable :: os(:)
2020-06-03 12:06:16 +02:00
! Memory allocation
2020-10-05 16:58:19 +02:00
maxS = min(nS,maxS)
allocate(X(nS),Y(nS),os(maxS))
2020-10-04 14:22:38 +02:00
! Compute oscillator strengths
2020-09-30 09:59:18 +02:00
2020-10-04 14:22:38 +02:00
os(:) = 0d0
2020-10-05 16:58:19 +02:00
if(spin_allowed) call oscillator_strength(nBas,nC,nO,nV,nR,nS,maxS,dipole_int,Omega,XpY,XmY,os)
2020-09-30 09:59:18 +02:00
! Print details about excitations
2020-06-03 12:06:16 +02:00
2020-10-05 16:58:19 +02:00
do ia=1,maxS
2020-06-03 12:06:16 +02:00
2020-06-05 22:34:32 +02:00
X(:) = 0.5d0*(XpY(ia,:) + XmY(ia,:))
Y(:) = 0.5d0*(XpY(ia,:) - XmY(ia,:))
2020-06-03 12:06:16 +02:00
2020-10-05 16:58:19 +02:00
! <S**2> values
if(spin_allowed) then
S2 = 0d0
else
S2 = 2d0
end if
2020-10-04 14:22:38 +02:00
print*,'-------------------------------------------------------------'
write(*,'(A15,I3,A2,F10.6,A3,A6,F6.4,A11,F6.4)') &
' Excitation n. ',ia,': ',Omega(ia)*HaToeV,' eV',' f = ',os(ia),' <S**2> = ',S2
print*,'-------------------------------------------------------------'
2020-06-03 12:06:16 +02:00
jb = 0
do j=nC+1,nO
do b=nO+1,nBas-nR
jb = jb + 1
2020-06-05 22:34:32 +02:00
if(abs(X(jb)) > thres_vec) write(*,'(I3,A4,I3,A3,F10.6)') j,' -> ',b,' = ',X(jb)/sqrt(2d0)
2020-06-03 12:06:16 +02:00
end do
end do
jb = 0
do j=nC+1,nO
do b=nO+1,nBas-nR
jb = jb + 1
2020-06-05 22:34:32 +02:00
if(abs(Y(jb)) > thres_vec) write(*,'(I3,A4,I3,A3,F10.6)') j,' <- ',b,' = ',Y(jb)/sqrt(2d0)
2020-06-03 12:06:16 +02:00
end do
end do
write(*,*)
end do
2020-10-03 22:16:38 +02:00
! Thomas-Reiche-Kuhn sum rule
write(*,'(A30,F10.6)') 'Thomas-Reiche-Kuhn sum rule = ',sum(os(:))
write(*,*)
2020-06-03 12:06:16 +02:00
end subroutine print_transition_vectors