The Elk Code
modfxcifc.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 module modfxcifc
7 
8 use libxcifc
9 
10 contains
11 
12 subroutine fxcifc(fxctype,n,rho,rhoup,rhodn,fxc,fxcuu,fxcud,fxcdd)
13 implicit none
14 ! mandatory arguments
15 integer, intent(in) :: fxctype(3),n
16 ! optional arguments
17 real(8), optional, intent(in) :: rho(n),rhoup(n),rhodn(n)
18 real(8), optional, intent(out) :: fxc(n),fxcuu(n),fxcud(n),fxcdd(n)
19 ! allocatable arrays
20 real(8), allocatable :: ra(:,:)
21 if (n < 1) then
22  write(*,*)
23  write(*,'("Error(fxcifc): n < 1 : ",I8)') n
24  write(*,*)
25  stop
26 end if
27 select case(abs(fxctype(1)))
28 case(0,1)
29 ! f_xc = 0
30  if (present(fxcuu).and.present(fxcud).and.present(fxcdd)) then
31  fxcuu(1:n)=0.d0
32  fxcud(1:n)=0.d0
33  fxcdd(1:n)=0.d0
34  else if (present(fxc)) then
35  fxc(1:n)=0.d0
36  else
37  goto 10
38  end if
39 case(3)
40 ! Perdew-Wang-Ceperley-Alder
41  if (present(rhoup).and.present(rhodn).and.present(fxcuu).and.present(fxcud) &
42  .and.present(fxcdd)) then
43 ! spin-polarised density
44  call fxc_pwca(n,rhoup,rhodn,fxcuu,fxcud,fxcdd)
45  else if (present(rho).and.present(fxc)) then
46 ! divide spin-unpolarised density into up and down
47  allocate(ra(n,4))
48  ra(1:n,1)=0.5d0*rho(1:n)
49  call fxc_pwca(n,ra(:,1),ra(:,1),ra(:,2),ra(:,3),ra(:,4))
50  fxc(1:n)=0.5d0*(ra(1:n,2)+ra(1:n,3))
51  deallocate(ra)
52  else
53  goto 10
54  end if
55 case(100)
56 ! libxc library functionals
57  if (present(rhoup).and.present(rhodn).and.present(fxcuu).and.present(fxcud) &
58  .and.present(fxcdd)) then
59 ! LSDA
60  call fxcifc_libxc(fxctype,n,rhoup=rhoup,rhodn=rhodn,fxcuu=fxcuu, &
61  fxcud=fxcud,fxcdd=fxcdd)
62  else if (present(rho).and.present(fxc)) then
63 ! LDA
64  call fxcifc_libxc(fxctype,n,rho=rho,fxc=fxc)
65  else
66  goto 10
67  end if
68 case default
69  write(*,*)
70  write(*,'("Error(fxcifc): response function unavailable for fxctype ",3I8)') &
71  fxctype
72  write(*,*)
73  stop
74 end select
75 return
76 10 continue
77 write(*,*)
78 write(*,'("Error(fxcifc): missing arguments for exchange-correlation type ",&
79  &3I6)') fxctype(:)
80 write(*,*)
81 stop
82 end subroutine
83 
84 end module
85 
subroutine fxcifc_libxc(fxctype, n, rho, rhoup, rhodn, fxc, fxcuu, fxcud, fxcdd)
subroutine fxc_pwca(n, rhoup, rhodn, fxcuu, fxcud, fxcdd)
Definition: fxc_pwca.f90:7
subroutine fxcifc(fxctype, n, rho, rhoup, rhodn, fxc, fxcuu, fxcud, fxcdd)
Definition: modfxcifc.f90:13