10
1
mirror of https://github.com/pfloos/quack synced 2024-12-23 21:03:47 +01:00
QuAcK/src/CI/UCIS.f90

136 lines
3.9 KiB
Fortran
Raw Normal View History

2020-10-05 23:00:56 +02:00
subroutine UCIS(spin_conserved,spin_flip,nBas,nC,nO,nV,nR,nS,ERI_aaaa,ERI_aabb,ERI_bbbb, &
2020-10-05 16:58:19 +02:00
dipole_int_aa,dipole_int_bb,eHF,cHF,S)
2020-09-24 14:39:37 +02:00
! Perform configuration interaction single calculation`
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: spin_conserved
logical,intent(in) :: spin_flip
integer,intent(in) :: nBas
integer,intent(in) :: nC(nspin)
integer,intent(in) :: nO(nspin)
integer,intent(in) :: nV(nspin)
integer,intent(in) :: nR(nspin)
integer,intent(in) :: nS(nspin)
double precision,intent(in) :: eHF(nBas,nspin)
2020-10-05 16:58:19 +02:00
double precision,intent(in) :: cHF(nBas,nBas,nspin)
double precision,intent(in) :: S(nBas,nBas)
2020-09-24 14:39:37 +02:00
double precision,intent(in) :: ERI_aaaa(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_aabb(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_bbbb(nBas,nBas,nBas,nBas)
2020-09-30 15:02:48 +02:00
double precision,intent(in) :: dipole_int_aa(nBas,nBas,ncart,nspin)
double precision,intent(in) :: dipole_int_bb(nBas,nBas,ncart,nspin)
2020-09-24 14:39:37 +02:00
! Local variables
logical :: dump_matrix = .false.
2020-10-05 16:58:19 +02:00
logical :: dump_trans = .false.
2020-09-24 14:39:37 +02:00
integer :: ispin
double precision :: lambda
integer :: nS_aa,nS_bb,nS_sc
double precision,allocatable :: A_sc(:,:)
2023-07-28 14:14:35 +02:00
double precision,allocatable :: Om_sc(:)
2020-09-24 14:39:37 +02:00
integer :: nS_ab,nS_ba,nS_sf
double precision,allocatable :: A_sf(:,:)
2023-07-28 14:14:35 +02:00
double precision,allocatable :: Om_sf(:)
2020-09-24 14:39:37 +02:00
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'| Configuration Interaction Singles |'
write(*,*)'************************************************'
write(*,*)
! Adiabatic connection scaling
lambda = 1d0
!----------------------------!
! Spin-conserved transitions !
!----------------------------!
if(spin_conserved) then
ispin = 1
! Memory allocation
nS_aa = nS(1)
nS_bb = nS(2)
nS_sc = nS_aa + nS_bb
2023-07-28 14:14:35 +02:00
allocate(A_sc(nS_sc,nS_sc),Om_sc(nS_sc))
2020-09-24 14:39:37 +02:00
2023-07-21 13:04:29 +02:00
call phULR_A(ispin,.false.,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,lambda,eHF,ERI_aaaa,ERI_aabb,ERI_bbbb,A_sc)
2020-09-24 14:39:37 +02:00
if(dump_matrix) then
print*,'CIS matrix (spin-conserved transitions)'
call matout(nS_sc,nS_sc,A_sc)
write(*,*)
endif
2023-07-28 14:14:35 +02:00
call diagonalize_matrix(nS_sc,A_sc,Om_sc)
2020-10-06 14:24:54 +02:00
A_sc(:,:) = transpose(A_sc)
2023-07-28 14:14:35 +02:00
call print_excitation_energies('UCIS',5,nS_sc,Om_sc)
call phULR_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nS_aa,nS_bb,nS_sc,dipole_int_aa,dipole_int_bb, &
cHF,S,Om_sc,A_sc,A_sc)
2020-09-24 14:39:37 +02:00
if(dump_trans) then
print*,'Spin-conserved CIS transition vectors'
call matout(nS_sc,nS_sc,A_sc)
write(*,*)
endif
2023-07-28 14:14:35 +02:00
deallocate(A_sc,Om_sc)
2020-09-24 14:39:37 +02:00
endif
!-----------------------!
! Spin-flip transitions !
!-----------------------!
if(spin_flip) then
ispin = 2
! Memory allocation
nS_ab = (nO(1) - nC(1))*(nV(2) - nR(2))
nS_ba = (nO(2) - nC(2))*(nV(1) - nR(1))
nS_sf = nS_ab + nS_ba
2023-07-28 14:14:35 +02:00
allocate(A_sf(nS_sf,nS_sf),Om_sf(nS_sf))
2020-09-30 15:02:48 +02:00
2023-07-21 13:04:29 +02:00
call phULR_A(ispin,.false.,nBas,nC,nO,nV,nR,nS_ab,nS_ba,nS_sf,lambda,eHF,ERI_aaaa,ERI_aabb,ERI_bbbb,A_sf)
2020-09-30 15:02:48 +02:00
2020-09-24 14:39:37 +02:00
if(dump_matrix) then
print*,'CIS matrix (spin-conserved transitions)'
call matout(nS_sf,nS_sf,A_sf)
write(*,*)
endif
2023-07-28 14:14:35 +02:00
call diagonalize_matrix(nS_sf,A_sf,Om_sf)
2020-10-06 14:24:54 +02:00
A_sf(:,:) = transpose(A_sf)
2023-07-28 14:14:35 +02:00
call print_excitation_energies('UCIS ',6,nS_sf,Om_sf)
call phULR_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nS_ab,nS_ba,nS_sf,dipole_int_aa,dipole_int_bb, &
cHF,S,Om_sf,A_sf,A_sf)
2020-09-24 14:39:37 +02:00
if(dump_trans) then
print*,'Spin-flip CIS transition vectors'
call matout(nS_sf,nS_sf,A_sf)
write(*,*)
endif
2023-07-28 14:14:35 +02:00
deallocate(A_sf,Om_sf)
2020-09-24 14:39:37 +02:00
endif
end subroutine