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 6module modfxcifc 7 8use libxcifc 9 10contains 11 12subroutine fxcifc(fxctype,n,rho,rhoup,rhodn,fxc,fxcuu,fxcud,fxcdd) 13implicit none 14! mandatory arguments 15integer, intent(in) :: fxctype(3),n 16! optional arguments 17real(8), optional, intent(in) :: rho(n),rhoup(n),rhodn(n) 18real(8), optional, intent(out) :: fxc(n),fxcuu(n),fxcud(n),fxcdd(n) 19! allocatable arrays 20real(8), allocatable :: ra(:,:) 21if (n.le.0) then 22 write(*,*) 23 write(*,'("Error(fxcifc): n <= 0 : ",I8)') n 24 write(*,*) 25 stop 26end if 27select case(abs(fxctype(1))) 28case(0,1) 29! f_xc = 0 30 if (present(fxcuu).and.present(fxcud).and.present(fxcdd)) then 31 fxcuu(:)=0.d0 32 fxcud(:)=0.d0 33 fxcdd(:)=0.d0 34 else if (present(fxc)) then 35 fxc(:)=0.d0 36 else 37 goto 10 38 end if 39case(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)=0.5d0*rho(:) 49 call fxc_pwca(n,ra(:,1),ra(:,1),ra(:,2),ra(:,3),ra(:,4)) 50 fxc(:)=0.5d0*(ra(:,2)+ra(:,3)) 51 deallocate(ra) 52 else 53 goto 10 54 end if 55case(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 68case default 69 write(*,*) 70 write(*,'("Error(fxcifc): response function unavailable for fxctype ",3I8)') & 71 fxctype 72 write(*,*) 73 stop 74end select 75return 7610 continue 77write(*,*) 78write(*,'("Error(fxcifc): missing arguments for exchange-correlation type ",& 79 &3I6)') fxctype(:) 80write(*,*) 81stop 82end subroutine 83 84end module 85 86