improved version of pw->wannier90 interface ...

git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2642 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
degironc 2005-12-22 14:01:00 +00:00
parent fe36de1c1b
commit b905dd4f58
1 changed files with 158 additions and 118 deletions

View File

@ -13,7 +13,7 @@ module wannier
integer, allocatable :: kpb(:,:) ! k+b (ik,ib) integer, allocatable :: kpb(:,:) ! k+b (ik,ib)
integer, allocatable :: g_kpb(:,:,:) ! G_k+b (ipol,ik,ib) integer, allocatable :: g_kpb(:,:,:) ! G_k+b (ipol,ik,ib)
integer, allocatable :: ig_(:,:) ! G_k+b (ipol,ik,ib) integer, allocatable :: ig_(:,:) ! G_k+b (ipol,ik,ib)
integer :: iun_nnkp, iun_mmn, iun_amn, nnbx, n_wannier integer :: iun_nnkp, iun_mmn, iun_amn, iun_band, nnbx, n_wannier
complex(DP), allocatable :: gf(:,:) ! guding_function(npwx,n_wannier) complex(DP), allocatable :: gf(:,:) ! guding_function(npwx,n_wannier)
real(DP), allocatable :: rw(:,:) ! rw(3,n_wannier) real(DP), allocatable :: rw(:,:) ! rw(3,n_wannier)
real(DP), allocatable :: alpha_w(:) ! alpha_w(n_wannier) real(DP), allocatable :: alpha_w(:) ! alpha_w(n_wannier)
@ -21,118 +21,124 @@ end module wannier
! !
!----------------------------------------------------------------------- !-----------------------------------------------------------------------
program pw2wannier90 program pw2wannier90
!----------------------------------------------------------------------- !-----------------------------------------------------------------------
! !
USE io_global, ONLY : stdout, ionode USE io_global, ONLY : stdout, ionode
USE mp_global, ONLY : mpime, kunit USE mp_global, ONLY : mpime, kunit
USE mp, ONLY : mp_bcast USE mp, ONLY : mp_bcast
USE cell_base, ONLY : at, bg USE cell_base, ONLY : at, bg
use lsda_mod, ONLY : nspin, isk use lsda_mod, ONLY : nspin, isk
use klist, ONLY : nkstot use klist, ONLY : nkstot
use ktetra, ONLY : k1, k2, k3, nk1, nk2, nk3 use ktetra, ONLY : k1, k2, k3, nk1, nk2, nk3
use io_files use io_files, ONLY : nd_nmbr, prefix, tmp_dir
! !
implicit none implicit none
integer :: ik, i, kunittmp integer :: ik, i, kunittmp
CHARACTER(LEN=4) :: spin_component CHARACTER(LEN=4) :: spin_component
integer :: ispinw CHARACTER(len=256) :: outdir
integer :: ispinw
namelist / inputpp / tmp_dir, prefix, spin_component
! namelist / inputpp / outdir, prefix, spin_component
call start_postproc (nd_nmbr) !
! call start_postproc (nd_nmbr)
! set default values for variables in namelist !
! ! set default values for variables in namelist
tmp_dir = './' !
prefix = ' ' outdir = './'
spin_component = 'none' prefix = ' '
! spin_component = 'none'
! reading the namelist inputpp !
! ! reading the namelist inputpp
read (5, inputpp) !
! read (5, inputpp)
! Check of namelist variables !
! ! Check of namelist variables
SELECT CASE ( TRIM( spin_component ) ) !
CASE ( 'up' ) tmp_dir = TRIM(outdir)
ispinw = 1 SELECT CASE ( TRIM( spin_component ) )
CASE ( 'down' ) CASE ( 'up' )
ispinw = 2 ispinw = 1
CASE DEFAULT CASE ( 'down' )
ispinw = 0 ispinw = 2
END SELECT CASE DEFAULT
! ispinw = 0
! Now allocate space for pwscf variables, read and check them. END SELECT
! !
call read_file ! Now allocate space for pwscf variables, read and check them.
call openfil_pp !
! call read_file
call read_nnkp call openfil_pp
! !
call compute_mmn call read_nnkp
! !
call compute_amn call compute_mmn
! !
call stop_pp call compute_amn
stop !
call write_band
!
call stop_pp
stop
end program pw2wannier90 end program pw2wannier90
! !
!----------------------------------------------------------------------- !-----------------------------------------------------------------------
subroutine read_nnkp subroutine read_nnkp
!----------------------------------------------------------------------- !-----------------------------------------------------------------------
! !
use kinds, only: DP use kinds, only: DP
use constants, only: eps8 use constants, only: eps8
use klist, only: nkstot use klist, only: nkstot
use cell_base, only : at, bg use cell_base, only : at, bg
use gvect, only : g,gg use gvect, only : g,gg
use parser, only : find_free_unit use parser, only : find_free_unit
use wannier use wannier
implicit none implicit none
real(DP) :: g_(3), gg_ real(DP) :: g_(3), gg_
integer :: ik, ib, ig, ipol, idum integer :: ik, ib, ig, ipol, idum
iun_nnkp = find_free_unit() iun_nnkp = find_free_unit()
open (unit=iun_nnkp, file='wannier.nnkp',form='formatted') open (unit=iun_nnkp, file='wannier.nnkp',form='formatted')
allocate ( nnb(nkstot) )
! get nnbx ! get nnbx
nnbx=0 nnbx=0
do ik=1,nkstot do ik=1,nkstot
read (iun_nnkp,*) nnb(ik) read (iun_nnkp,*) nnb(ik)
nnbx = max (nnbx, nnb(ik) ) nnbx = max (nnbx, nnb(ik) )
do ib = 1, nnb(ik) do ib = 1, nnb(ik)
read(iun_nnkp,*) read(iun_nnkp,*)
end do end do
end do end do
! allocate ! allocate
allocate ( kpb(nkstot,nnbx), g_kpb(3,nkstot,nnbx),ig_(nkstot,nnbx) ) allocate ( kpb(nkstot,nnbx), g_kpb(3,nkstot,nnbx),ig_(nkstot,nnbx) )
! read ! read
rewind(iun_nnkp) rewind(iun_nnkp)
do ik=1,nkstot do ik=1,nkstot
read (iun_nnkp,*) idum read (iun_nnkp,*) idum
if (idum /= nnb(ik)) call errore('pw2w','ekkekka ',1) if (idum /= nnb(ik)) call errore('pw2w','ekkekka ',1)
do ib = 1, nnb(ik) do ib = 1, nnb(ik)
read(iun_nnkp,*) idum, kpb(ik,ib), (g_kpb(ipol,ik,ib), ipol =1,3) read(iun_nnkp,*) idum, kpb(ik,ib), (g_kpb(ipol,ik,ib), ipol =1,3)
g_(:) = REAL( g_kpb(:,ik,ib) ) g_(:) = REAL( g_kpb(:,ik,ib) )
call trnvecc (g_, bg, at, +1) call trnvect (g_, at, bg, 1)
gg_ = g_(1)*g_(1) + g_(2)*g_(2) + g_(3)*g_(3) gg_ = g_(1)*g_(1) + g_(2)*g_(2) + g_(3)*g_(3)
ig_(ik,ib) = 0 ig_(ik,ib) = 0
ig = 1 ig = 1
do while (gg(ig) < gg_) do while (gg(ig) <= gg_ + eps8)
if ( (abs(g(1,ig)-g_(1)) < eps8) .and. & if ( (abs(g(1,ig)-g_(1)) < eps8) .and. &
(abs(g(2,ig)-g_(2)) < eps8) .and. & (abs(g(2,ig)-g_(2)) < eps8) .and. &
(abs(g(3,ig)-g_(3)) < eps8) ) ig_(ik,ib) = ig (abs(g(3,ig)-g_(3)) < eps8) ) ig_(ik,ib) = ig
ig= ig +1 ig= ig +1
end do end do
end do ! write (*,'(4i6,3f10.6)') ik, ib, kpb(ik,ib), ig_(ik,ib), g_(:)
end do end do
end do
close (iun_nnkp)
close (iun_nnkp) return
return
end subroutine read_nnkp end subroutine read_nnkp
! !
@ -154,6 +160,7 @@ subroutine compute_mmn
complex(DP), allocatable :: phase(:), aux(:), evcq(:,:) complex(DP), allocatable :: phase(:), aux(:), evcq(:,:)
integer, allocatable :: igkq(:) integer, allocatable :: igkq(:)
complex(DP) :: mmn, ZDOTC complex(DP) :: mmn, ZDOTC
real(DP)::aa
allocate( phase(nrxxs), aux(npwx), evcq(npwx,nbnd), igkq(npwx) ) allocate( phase(nrxxs), aux(npwx), evcq(npwx,nbnd), igkq(npwx) )
@ -165,8 +172,10 @@ subroutine compute_mmn
mmn_tot = mmn_tot + nnb(ik) * nbnd * nbnd mmn_tot = mmn_tot + nnb(ik) * nbnd * nbnd
end do end do
write (*,*) "MMN"
write (iun_mmn,*) mmn_tot write (iun_mmn,*) mmn_tot
do ik=1,nkstot do ik=1,nkstot
write (*,*) ik
! read wfc at k ! read wfc at k
call davcio (evc, nwordwfc, iunwfc, ik, -1 ) call davcio (evc, nwordwfc, iunwfc, ik, -1 )
call gk_sort (xk(1,ik), ngm, g, ecutwfc / tpiba2, npw, igk, g2kin) call gk_sort (xk(1,ik), ngm, g, ecutwfc / tpiba2, npw, igk, g2kin)
@ -175,7 +184,7 @@ subroutine compute_mmn
ikp = kpb(ik,ib) ikp = kpb(ik,ib)
! read wfc at k+b ! read wfc at k+b
call davcio (evcq, nwordwfc, iunwfc, ikp, -1 ) call davcio (evcq, nwordwfc, iunwfc, ikp, -1 )
call gk_sort (xk(1,ikp), ngm, g, ecutwfc / tpiba2, npw, igkq, g2kin) call gk_sort (xk(1,ikp), ngm, g, ecutwfc / tpiba2, npwq, igkq, g2kin)
! compute the phase ! compute the phase
phase(:) = (0.d0,0.d0) phase(:) = (0.d0,0.d0)
if ( ig_(ik,ib)>0) phase( nls(ig_(ik,ib)) ) = (1.d0,0.d0) if ( ig_(ik,ib)>0) phase( nls(ig_(ik,ib)) ) = (1.d0,0.d0)
@ -186,15 +195,18 @@ subroutine compute_mmn
psic(nls (igk (1:npw) ) ) = evc (1:npw, ibnd) psic(nls (igk (1:npw) ) ) = evc (1:npw, ibnd)
call cft3s (psic, nr1s, nr2s, nr3s, nrx1s, nrx2s, nrx3s, +2) call cft3s (psic, nr1s, nr2s, nr3s, nrx1s, nrx2s, nrx3s, +2)
psic(1:nrxxs) = psic(1:nrxxs) * phase(1:nrxxs) psic(1:nrxxs) = psic(1:nrxxs) * phase(1:nrxxs)
call cft3s (psic, nr1s, nr2s, nr3s, nrx1s, nrx2s, nrx3s, +2) call cft3s (psic, nr1s, nr2s, nr3s, nrx1s, nrx2s, nrx3s, -2)
aux(1:npwq) = psic(nls (igkq(1:npwq) ) ) aux(1:npwq) = psic(nls (igkq(1:npwq) ) )
aa = 0.d0
do jbnd=1,nbnd do jbnd=1,nbnd
mmn = ZDOTC ( npwq, aux,1,evcq(1,jbnd),1) mmn = ZDOTC ( npwq, aux,1,evcq(1,jbnd),1)
call reduce(mmn,2) call reduce(mmn,2)
write (iun_mmn,'(7i6)') ibnd, jbnd, ik, ikp, & write (iun_mmn,'(7i5)') ibnd, jbnd, ik, ikp, &
(g_kpb(ipol,ik,ib), ipol=1,3) (g_kpb(ipol,ik,ib), ipol=1,3)
write (iun_mmn,*) mmn write (iun_mmn,'(2f18.12)') mmn
aa = aa + abs(mmn)**2
end do end do
! write (*,*) ik,ib,ibnd, aa
end do end do
end do end do
end do end do
@ -222,20 +234,22 @@ subroutine compute_amn
call read_gf_definition call read_gf_definition
iun_amn = find_free_unit() iun_amn = find_free_unit()
open (unit=iun_amn, file='wannier.mmn',form='formatted') open (unit=iun_amn, file='wannier.amn',form='formatted')
amn_tot = nkstot * nbnd * n_wannier amn_tot = nkstot * nbnd * n_wannier
write (*,*) "AMN"
write (iun_amn,*) amn_tot write (iun_amn,*) amn_tot
do ik=1,nkstot do ik=1,nkstot
write (*,*) ik
call davcio (evc, nwordwfc, iunwfc, ik, -1 ) call davcio (evc, nwordwfc, iunwfc, ik, -1 )
call gk_sort (xk(1,ik), ngm, g, ecutwfc / tpiba2, npw, igk, g2kin) call gk_sort (xk(1,ik), ngm, g, ecutwfc / tpiba2, npw, igk, g2kin)
call generate_guiding_functions call generate_guiding_functions(ik)
do ibnd=1,nbnd do ibnd=1,nbnd
do iw=1,n_wannier do iw=1,n_wannier
amn = ZDOTC(npw,evc(1,ibnd),1,gf(1,iw),1) amn = ZDOTC(npw,evc(1,ibnd),1,gf(1,iw),1)
call reduce(amn,2) call reduce(amn,2)
write(iun_amn,*) ibnd,iw,ik, amn write(iun_amn,'(3i5,2f18.12)') ibnd,iw,ik, amn
end do end do
end do end do
end do end do
@ -257,24 +271,50 @@ subroutine read_gf_definition
do iw=1,n_wannier do iw=1,n_wannier
read (*,*) (rw(i,iw),i=1,3), alpha_w(iw) read (*,*) (rw(i,iw),i=1,3), alpha_w(iw)
end do end do
! convert them in carthesian coordinated in unit of alat
! rw(:,:) = rw(:,:) / celldm(1)/0.529177d0
return return
end subroutine read_gf_definition end subroutine read_gf_definition
! !
subroutine generate_guiding_functions subroutine generate_guiding_functions(ik)
use constants, only : tpi use constants, only : pi, tpi, fpi
use wvfct, only : npw, g2kin use wvfct, only : npw, g2kin
use gvect, only : g use gvect, only : g
use cell_base, ONLY :tpiba2, omega
use wannier use wannier
implicit none implicit none
integer iw, ig integer iw, ig, ik
real(DP) :: arg real(DP) :: arg, anorm, fac, alpha_w2
complex(DP) :: ZDOTC
do iw =1, n_wannier do iw =1, n_wannier
fac = ( sqrt(fpi)*alpha_w(iw) )**(1.5) / sqrt(omega)
alpha_w2 = alpha_w(iw)**2
do ig=1,npw do ig=1,npw
arg = tpi * (g(1,ig)*rw(1,iw)+g(2,ig)*rw(2,iw)+g(3,ig)*rw(3,iw)) arg = tpi * (g(1,ig)*rw(1,iw)+g(2,ig)*rw(2,iw)+g(3,ig)*rw(3,iw))
gf(ig,iw) = exp(-alpha_w(iw)*g2kin(ig)) * & gf(ig,iw) = fac * exp(-0.5d0*alpha_w2*g2kin(ig)*tpiba2) * &
CMPLX(cos(arg),sin(arg)) CMPLX(cos(arg),sin(arg))
end do end do
anorm = ZDOTC(npw,gf(1,iw),1,gf(1,iw),1)
write (*,*) iw, ik, anorm
end do end do
return return
end subroutine generate_guiding_functions end subroutine generate_guiding_functions
subroutine write_band
use wvfct, only : nbnd, et
use klist, only : nkstot
use constants, only: rytoev
use parser, only : find_free_unit
use wannier
implicit none
integer ik, ibnd
iun_band = find_free_unit()
open (unit=iun_band, file='BAND.dat',form='formatted')
do ik=1,nkstot
do ibnd=1,nbnd
write (iun_band,'(2i5,f18.12)') ibnd, ik, et(ibnd,ik)*rytoev
end do
end do
return
end subroutine write_band