2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! Copyright (C) 2001 PWSCF group
|
|
|
|
! This file is distributed under the terms of the
|
|
|
|
! GNU General Public License. See the file `License'
|
|
|
|
! in the root directory of the present distribution,
|
|
|
|
! or http://www.gnu.org/copyleft/gpl.txt .
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine addusldos (ldos, becsum1)
|
2003-01-20 05:58:50 +08:00
|
|
|
!----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine adds to the local DOS the part which is due to
|
|
|
|
! the US augmentation.
|
|
|
|
!
|
|
|
|
#include "machine.h"
|
2003-02-08 00:04:36 +08:00
|
|
|
use pwcom
|
2003-01-20 05:58:50 +08:00
|
|
|
implicit none
|
2003-02-08 00:04:36 +08:00
|
|
|
complex(kind=DP) :: ldos (nrxx, nspin)
|
2003-01-20 05:58:50 +08:00
|
|
|
! local density of states
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
real(kind=DP) :: becsum1 ( (nhm * (nhm + 1) ) / 2, nat, nspin)
|
2003-01-20 05:58:50 +08:00
|
|
|
! input: the becsum1 ter
|
|
|
|
!
|
|
|
|
! here the local variables
|
|
|
|
!
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
integer :: ig, na, nt, ih, jh, ijh, is, ir
|
2003-01-20 05:58:50 +08:00
|
|
|
! counter on G vectors
|
|
|
|
! counter on atoms
|
|
|
|
! counter on atomic types
|
|
|
|
! counter on beta functions
|
|
|
|
! counter on beta functions
|
|
|
|
! composite index ih jh
|
|
|
|
! counter on spin
|
|
|
|
! counter on mesh points
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
real(kind=DP),pointer :: ylmk0 (:,:), qmod (:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! the spherical harmonics
|
|
|
|
! the modulus of G
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
complex(kind=DP), allocatable :: qg (:), aux (:,:)
|
2003-01-20 05:58:50 +08:00
|
|
|
! auxiliary variable for FFT
|
|
|
|
! auxiliary variable for rho(G)
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (aux ( ngm , nspin))
|
|
|
|
allocate (qg ( nrxx))
|
|
|
|
allocate (qmod( ngm))
|
|
|
|
allocate (ylmk0 ( ngm , lqx * lqx))
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call setv (2 * ngm * nspin, 0.d0, aux, 1)
|
|
|
|
call ylmr2 (lqx * lqx, ngm, g, gg, ylmk0)
|
|
|
|
do ig = 1, ngm
|
|
|
|
qmod (ig) = sqrt (gg (ig) )
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
do nt = 1, ntyp
|
|
|
|
if (tvanp (nt) ) then
|
|
|
|
ijh = 0
|
|
|
|
do ih = 1, nh (nt)
|
|
|
|
do jh = ih, nh (nt)
|
|
|
|
call qvan2 (ngm, ih, jh, nt, qmod, qgm, ylmk0)
|
|
|
|
ijh = ijh + 1
|
|
|
|
do na = 1, nat
|
|
|
|
if (ityp (na) .eq.nt) then
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! Multiply becsum and qg with the correct structure factor
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
do is = 1, nspin
|
|
|
|
do ig = 1, ngm
|
2003-01-20 05:58:50 +08:00
|
|
|
aux (ig, is) = aux (ig, is) + qgm (ig) * becsum1 (ijh, na, &
|
|
|
|
is) * (eigts1 (ig1 (ig), na) * eigts2 (ig2 (ig), na) &
|
|
|
|
* eigts3 (ig3 (ig), na) )
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! convert aux to real space and adds to the charge density
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
if (okvan) then
|
|
|
|
do is = 1, nspin
|
|
|
|
call setv (2 * nrxx, 0.d0, qg, 1)
|
|
|
|
do ig = 1, ngm
|
|
|
|
qg (nl (ig) ) = aux (ig, is)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
enddo
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
call cft3 (qg, nr1, nr2, nr3, nrx1, nrx2, nrx3, 1)
|
|
|
|
do ir = 1, nrxx
|
|
|
|
ldos (ir, is) = ldos (ir, is) + DREAL (qg (ir) )
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
|
|
|
|
enddo
|
|
|
|
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
deallocate (ylmk0)
|
|
|
|
deallocate (qmod)
|
|
|
|
deallocate (qg)
|
|
|
|
deallocate (aux)
|
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine addusldos
|