4
1
mirror of https://github.com/pfloos/quack synced 2024-06-02 11:25:28 +02:00
quack/src/CI/RCIS.f90

121 lines
2.9 KiB
Fortran
Raw Normal View History

subroutine RCIS(dotest,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
logical,intent(in) :: dotest
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)
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.
logical :: dRPA = .false.
2019-03-19 10:13:33 +01:00
integer :: ispin
integer :: maxS = 10
2020-01-08 10:17:19 +01:00
double precision :: lambda
double precision,allocatable :: A(:,:),Om(:)
2019-03-19 10:13:33 +01:00
! Hello world
write(*,*)
2023-11-22 14:54:05 +01:00
write(*,*)'******************************'
write(*,*)'* Restricted CIS Calculation *'
write(*,*)'******************************'
2019-03-19 10:13:33 +01:00
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),Om(nS))
2019-03-19 10:13:33 +01:00
! Compute CIS matrix
if(singlet) then
2019-03-19 10:13:33 +01:00
ispin = 1
call phLR_A(ispin,dRPA,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(*,*)
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
call diagonalize_matrix(nS,A,Om)
2023-11-22 10:07:23 +01:00
call print_excitation_energies('CIS@RHF','singlet',nS,Om)
2023-07-28 14:14:35 +02:00
call phLR_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,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(*,*)
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
! Compute CIS(D) correction
maxS = min(maxS,nS)
2023-07-28 14:14:35 +02:00
if(doCIS_D) call CIS_D(ispin,nBas,nC,nO,nV,nR,nS,maxS,eHF,ERI,Om(1:maxS),A(:,1:maxS))
! Testing zone
if(dotest) then
2023-11-14 14:31:27 +01:00
call dump_test_value('R','CIS singlet excitation energy',Om(1))
end if
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
if(triplet) then
2019-03-19 10:13:33 +01:00
ispin = 2
call phLR_A(ispin,dRPA,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(*,*)
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
call diagonalize_matrix(nS,A,Om)
2023-11-22 10:07:23 +01:00
call print_excitation_energies('CIS@RHF','triplet',nS,Om)
2023-07-28 14:14:35 +02:00
call phLR_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,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(*,*)
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
! Compute CIS(D) correction
maxS = min(maxS,nS)
2023-07-28 14:14:35 +02:00
if(doCIS_D) call CIS_D(ispin,nBas,nC,nO,nV,nR,nS,maxS,eHF,ERI,Om(1:maxS),A(:,1:maxS))
! Testing zone
if(dotest) then
2023-11-14 14:31:27 +01:00
call dump_test_value('R','CIS triplet excitation energy',Om(1))
end if
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
2023-07-12 22:37:04 +02:00
end subroutine