4
1
mirror of https://github.com/pfloos/quack synced 2024-11-07 22:53:59 +01:00
quack/src/QuAcK/CIS.f90

88 lines
2.0 KiB
Fortran
Raw Normal View History

2019-03-19 10:13:33 +01:00
subroutine CIS(singlet_manifold,triplet_manifold, &
nBas,nC,nO,nV,nR,nS,ERI,eHF)
! Perform configuration interaction single calculation`
implicit none
include 'parameters.h'
! Input variables
2020-01-15 22:29:43 +01:00
logical,intent(in) :: singlet_manifold
logical,intent(in) :: triplet_manifold
2019-03-19 10:13:33 +01:00
integer,intent(in) :: nBas,nC,nO,nV,nR,nS
2020-01-15 22:29:43 +01:00
double precision,intent(in) :: eHF(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
2019-03-19 10:13:33 +01:00
! Local variables
logical :: dump_matrix = .false.
logical :: dump_trans = .false.
integer :: ispin
2020-01-08 10:17:19 +01:00
double precision :: lambda
2019-03-19 10:13:33 +01:00
double precision,allocatable :: A(:,:),Omega(:)
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'| Configuration Interaction Singles |'
write(*,*)'************************************************'
write(*,*)
2020-01-08 10:17:19 +01:00
! Adiabatic connection scaling
lambda = 1d0
2019-03-19 10:13:33 +01:00
! Memory allocation
allocate(A(nS,nS),Omega(nS))
! Compute CIS matrix
if(singlet_manifold) then
ispin = 1
2020-01-15 22:29:43 +01:00
call linear_response_A_matrix(ispin,.false.,nBas,nC,nO,nV,nR,nS,lambda,eHF,ERI,A)
2019-03-19 10:13:33 +01:00
if(dump_matrix) then
print*,'CIS matrix (singlet state)'
call matout(nS,nS,A)
write(*,*)
endif
call diagonalize_matrix(nS,A,Omega)
2020-04-26 16:07:45 +02:00
call print_excitation('CIS ',ispin,nS,Omega)
2019-03-19 10:13:33 +01:00
if(dump_trans) then
print*,'Singlet CIS transition vectors'
call matout(nS,nS,A)
write(*,*)
endif
endif
if(triplet_manifold) then
ispin = 2
2020-01-15 22:29:43 +01:00
call linear_response_A_matrix(ispin,.false.,nBas,nC,nO,nV,nR,nS,lambda,eHF,ERI,A)
2019-03-19 10:13:33 +01:00
if(dump_matrix) then
print*,'CIS matrix (triplet state)'
call matout(nS,nS,A)
write(*,*)
endif
call diagonalize_matrix(nS,A,Omega)
2020-04-26 16:07:45 +02:00
call print_excitation('CIS ',ispin,nS,Omega)
2019-03-19 10:13:33 +01:00
if(dump_trans) then
print*,'Triplet CIS transition vectors'
call matout(nS,nS,A)
write(*,*)
endif
endif
end subroutine CIS