TABLE OF CONTENTS
ABINIT/upfheader2abi [ Functions ]
NAME
upfheader2abi
FUNCTION
This routine wraps a call to a PWSCF module, which reads in a UPF (PWSCF / Espresso) format pseudopotential, then transfers data for the HEADER of abinit psps only!
COPYRIGHT
Copyright (C) 2009-2018 ABINIT group (MJV) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
filpsp = name of file with UPF data
OUTPUT
pspxc = index of xc functional for this pseudo lmax_ = maximal angular momentum znucl = charge of species nucleus zion = valence charge n1xccc = default number of points. Set to 0 if no nlcc is present nproj_l= number of projectors for each channel nprojso_l= number of projectors for each channel for SO correction projectors
SIDE EFFECTS
NOTES
PARENTS
inpspheads
CHILDREN
set_dft_from_indices,set_dft_from_name
SOURCE
42 #if defined HAVE_CONFIG_H 43 #include "config.h" 44 #endif 45 46 #include "abi_common.h" 47 48 49 subroutine upfheader2abi (filpsp, znucl, zion, pspxc, lmax_, n1xccc, nproj_l, nprojso_l) 50 51 use defs_basis 52 use m_profiling_abi 53 use m_errors 54 use m_atomdata 55 use pseudo_pwscf ! pwscf module with all data explicit! 56 57 use m_io_tools, only : open_file 58 59 !This section has been created automatically by the script Abilint (TD). 60 !Do not modify the following lines by hand. 61 #undef ABI_FUNC 62 #define ABI_FUNC 'upfheader2abi' 63 use interfaces_11_qespresso_ext 64 use interfaces_57_iopsp_parser, except_this_one => upfheader2abi 65 !End of the abilint section 66 67 implicit none 68 69 !Arguments ------------------------------- 70 character(len=fnlen), intent(in) :: filpsp 71 integer, intent(inout) :: n1xccc 72 integer, intent(out) :: pspxc, lmax_ 73 real(dp), intent(out) :: znucl, zion 74 !arrays 75 integer, intent(out) :: nproj_l(0:3) 76 integer, intent(out) :: nprojso_l(1:3) 77 78 !Local variables ------------------------- 79 integer :: iproj, ll, iunit 80 character(len=500) :: msg 81 type(atomdata_t) :: atom 82 83 ! ********************************************************************* 84 85 !call pwscf routine for reading in UPF 86 if (open_file(filpsp, msg, newunit=iunit, status='old',form='formatted') /= 0) then 87 MSG_ERROR(msg) 88 end if 89 90 !read in psp data to static data in pseudo module, for ipsx == 1 91 call read_pseudo(1,iunit) 92 close (iunit) 93 94 !copy over to abinit internal arrays and vars 95 call upfxc2abi(dft(1), pspxc) 96 lmax_ = lmax(1) 97 call atomdata_from_symbol(atom,psd(1)) 98 znucl = atom%znucl 99 zion = zp(1) 100 101 nproj_l = 0 102 do iproj = 1, nbeta(1) 103 ll = lll(iproj,1) 104 nproj_l(ll) = nproj_l(ll) + 1 105 end do 106 107 nprojso_l = 0 !FIXME deal with so 108 !do iproj = 1, nbeta(1) 109 !nprojso_l(ll+1) = nprojso_l(ll+1) + 1 110 !end do 111 112 if (.not. nlcc(1)) n1xccc = 0 113 114 end subroutine upfheader2abi
ABINIT/upfxc2abi [ Functions ]
NAME
upfxc2abi
FUNCTION
This routine wraps a call to an OCTOPUS module, which reformats a UPF (PWSCF / Espresso) string describing XC functionals, and returns the abinit internal code pspxc
COPYRIGHT
Copyright (C) 2009-2018 ABINIT group (MJV) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
dft = string with x/c functionals from PWSCF format
OUTPUT
pspxc = index of xc functional for this pseudo
SIDE EFFECTS
NOTES
FIXME: extend to more functionals with libxc Could be included in separate module, eg read_upf_pwscf or funct_pwscf Left without defs_basis or calls to abinit routines ON PURPOSE
PARENTS
upf2abinit,upfheader2abi
CHILDREN
set_dft_from_indices,set_dft_from_name
SOURCE
153 subroutine upfxc2abi(dft, pspxc) 154 155 use defs_basis 156 use funct_pwscf ! pwscf module for naming xc functionals 157 use m_profiling_abi 158 use m_errors 159 160 !This section has been created automatically by the script Abilint (TD). 161 !Do not modify the following lines by hand. 162 #undef ABI_FUNC 163 #define ABI_FUNC 'upfxc2abi' 164 !End of the abilint section 165 166 implicit none 167 168 !Arguments ------------------------------- 169 character(len=20), intent(in) :: dft 170 integer, intent(out) :: pspxc 171 172 !Local variables ------------------------- 173 integer :: iexch,icorr,igcx,igcc 174 integer :: totalindex, offset 175 176 ! ********************************************************************* 177 !extract from char*20 :: dft(:) 178 !### The following has been copied from pwscf src/Modules/upf_to_internal.f90: 179 !workaround for rrkj format - it contains the indices, not the name 180 if ( dft(1:6)=='INDEX:') then 181 read( dft(7:10), '(4i1)') iexch,icorr,igcx,igcc 182 call set_dft_from_indices(iexch,icorr,igcx,igcc) 183 else 184 call set_dft_from_name( dft ) 185 iexch = get_iexch() 186 icorr = get_icorr() 187 igcx = get_igcx() 188 igcc = get_igcc() 189 end if 190 !reset dft string to avoid stray spaces 191 call set_dft_from_indices(iexch,icorr,igcx,igcc) 192 write(std_out,'(a)') ' upf2abinit: XC string from pseudopotential is :' 193 write(std_out,'(3a)') '>', dft, '<' 194 195 offset = 100 196 totalindex = offset*offset*offset*iexch + offset*offset*icorr + offset*igcx + igcc 197 select case (totalindex) 198 case (00000000) !(" NOX NOC NOGX NOGC") ! no xc 199 pspxc = 0 200 case (01010000) !(" SLA PZ NOGX NOGC") ! slater exchange + Perdew Zunger 201 pspxc = 2 202 case (01050000) !(" SLA WIG NOGX NOGC") ! slater exchange + Wigner corr 203 pspxc = 4 204 case (01060000) !(" SLA HL NOGX NOGC") ! Hedin + Lundqvist 205 pspxc = 5 206 case (02000000) !(" SL1 NOC NOGX NOGC") ! full slater exchange 207 pspxc = 6 208 case (01040000) !(" SLA PW NOGX NOGC") ! slater exchange + Perdew Wang 209 pspxc = 7 210 case (01000000) !(" SLA NOC NOGX NOGC") ! Perdew Wang + no corr 211 pspxc = 8 212 case (01040304) !(" SLA PW PBX PBC") ! LDA + PBE GGA 213 pspxc = 11 ! PBE 214 case (01000300) !(" SLA NOC PBX NOGC") ! exchange part of PBE GGA 215 pspxc = 12 216 case (01040404) !(" SLA PW RPB PBC") ! rev PBE 217 pspxc = 14 218 case (00000505) !(" NOX NOC HTCH HTCH") ! HTCH 120 219 pspxc = 17 220 case (01030103) !(" SLA LYP B88 BLYP") ! BLYP 221 pspxc = -106131 222 case (01040101) !(" SLA PW B88 P86") ! BP86 223 pspxc = -106132 224 case (00030603) !(" NOX LYP OPTX BLYP") ! OLYP 225 pspxc = -110131 226 ! FIXME: important cases left to be patched with libxc: 227 ! vosko wilkins nusair 228 ! ortiz ballone 229 ! pbe0 230 ! Gunnarson-Lunqvist 231 ! make general approach: check gradient parts first, then lda. 232 ! event. check if they are consistent. 233 case default 234 MSG_ERROR('upf2abinit: XC functional not recognized') 235 end select 236 237 end subroutine upfxc2abi