2020-10-01 23:05:23 +02:00
|
|
|
subroutine CIS(singlet,triplet,doCIS_D,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF)
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
! Perform configuration interaction single calculation`
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
include 'parameters.h'
|
|
|
|
|
|
|
|
! Input variables
|
|
|
|
|
2020-10-01 23:05:23 +02:00
|
|
|
logical,intent(in) :: singlet
|
|
|
|
logical,intent(in) :: triplet
|
|
|
|
logical,intent(in) :: doCIS_D
|
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)
|
2020-09-28 21:25:25 +02:00
|
|
|
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
! Local variables
|
|
|
|
|
|
|
|
logical :: dump_matrix = .false.
|
|
|
|
logical :: dump_trans = .false.
|
|
|
|
integer :: ispin
|
2020-10-01 23:05:23 +02:00
|
|
|
integer :: maxS = 10
|
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
|
|
|
|
|
2020-10-01 23:05:23 +02:00
|
|
|
if(singlet) then
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
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)
|
2020-10-03 22:16:38 +02:00
|
|
|
call print_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,Omega,transpose(A),transpose(A))
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
if(dump_trans) then
|
|
|
|
print*,'Singlet CIS transition vectors'
|
|
|
|
call matout(nS,nS,A)
|
|
|
|
write(*,*)
|
|
|
|
endif
|
|
|
|
|
2020-10-01 23:05:23 +02:00
|
|
|
! Compute CIS(D) correction
|
|
|
|
|
|
|
|
maxS = min(maxS,nS)
|
|
|
|
if(doCIS_D) call D_correction(ispin,nBas,nC,nO,nV,nR,nS,maxS,eHF,ERI,Omega(1:maxS),A(:,1:maxS))
|
|
|
|
|
2019-03-19 10:13:33 +01:00
|
|
|
endif
|
|
|
|
|
2020-10-01 23:05:23 +02:00
|
|
|
if(triplet) then
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
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)
|
2020-10-03 22:16:38 +02:00
|
|
|
call print_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int,Omega,transpose(A),transpose(A))
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
if(dump_trans) then
|
|
|
|
print*,'Triplet CIS transition vectors'
|
|
|
|
call matout(nS,nS,A)
|
|
|
|
write(*,*)
|
|
|
|
endif
|
|
|
|
|
2020-10-01 23:05:23 +02:00
|
|
|
! Compute CIS(D) correction
|
|
|
|
|
|
|
|
maxS = min(maxS,nS)
|
|
|
|
if(doCIS_D) call D_correction(ispin,nBas,nC,nO,nV,nR,nS,maxS,eHF,ERI,Omega(1:maxS),A(:,1:maxS))
|
|
|
|
|
2019-03-19 10:13:33 +01:00
|
|
|
endif
|
|
|
|
|
|
|
|
end subroutine CIS
|