4
1
mirror of https://github.com/pfloos/quack synced 2024-06-27 23:52:19 +02:00
quack/src/eDFT/RMFL20_lda_correlation_energy.f90

67 lines
1.7 KiB
Fortran
Raw Normal View History

2020-03-15 20:33:18 +01:00
subroutine RMFL20_lda_correlation_energy(nEns,wEns,nGrid,weight,rho,Ec)
2019-03-13 11:07:31 +01:00
2020-03-15 20:33:18 +01:00
! Compute the restricted version of the Marut-Fromager-Loos weight-dependent correlation functional
! The RMFL20 is a two-state, single-weight correlation functional for spin-unpolarized systems
2019-03-13 11:07:31 +01:00
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: nEns
double precision,intent(in) :: wEns(nEns)
integer,intent(in) :: nGrid
double precision,intent(in) :: weight(nGrid)
2020-03-17 11:50:05 +01:00
double precision,intent(in) :: rho(nGrid)
2019-03-13 11:07:31 +01:00
! Local variables
2020-03-29 18:34:09 +02:00
logical :: LDA_centered = .true.
2020-03-17 11:50:05 +01:00
integer :: iEns
double precision :: EcLDA
2020-03-15 20:33:18 +01:00
double precision,allocatable :: aMFL(:,:)
2020-03-17 11:50:05 +01:00
double precision,allocatable :: EceLDA(:)
2019-03-13 11:07:31 +01:00
! Output variables
2020-03-17 11:50:05 +01:00
double precision :: Ec
2019-03-13 11:07:31 +01:00
! Allocation
2020-03-17 11:50:05 +01:00
allocate(aMFL(3,nEns),EceLDA(nEns))
2019-03-13 11:07:31 +01:00
! Parameters for weight-dependent LDA correlation functional
2020-03-15 20:33:18 +01:00
aMFL(1,1) = -0.0238184d0
aMFL(2,1) = +0.00540994d0
aMFL(3,1) = +0.0830766d0
2019-03-13 11:07:31 +01:00
2020-03-15 20:33:18 +01:00
aMFL(1,2) = -0.0144633d0
aMFL(2,2) = -0.0506019d0
aMFL(3,2) = +0.0331417d0
2019-03-13 11:07:31 +01:00
2020-03-17 11:50:05 +01:00
! Compute correlation energy for ground and doubly-excited states
2019-03-13 11:07:31 +01:00
do iEns=1,nEns
2020-03-27 20:46:13 +01:00
call restricted_elda_correlation_energy(aMFL(:,iEns),nGrid,weight(:),rho(:),EceLDA(iEns))
2019-03-13 11:07:31 +01:00
end do
! LDA-centered functional
2020-03-17 11:50:05 +01:00
call RVWN5_lda_correlation_energy(nGrid,weight(:),rho(:),EcLDA)
2019-03-13 11:07:31 +01:00
2020-03-29 18:34:09 +02:00
if(LDA_centered) then
do iEns=1,nEns
EceLDA(iEns) = EceLDA(iEns) + EcLDA - EceLDA(1)
end do
end if
2019-03-13 11:07:31 +01:00
! Weight-denpendent functional for ensembles
2020-03-17 11:50:05 +01:00
Ec = 0d0
2019-03-13 11:07:31 +01:00
do iEns=1,nEns
2020-03-17 11:50:05 +01:00
Ec = Ec + wEns(iEns)*EceLDA(iEns)
2019-03-13 11:07:31 +01:00
end do
2020-03-15 20:33:18 +01:00
end subroutine RMFL20_lda_correlation_energy