2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-01-28 20:28:11 +08:00
|
|
|
! Copyright (C) 2003 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-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine cg_setup
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#include "machine.h"
|
2004-01-23 23:08:03 +08:00
|
|
|
USE kinds, only: DP
|
2003-01-20 05:58:50 +08:00
|
|
|
use pwcom
|
2004-04-28 18:25:36 +08:00
|
|
|
USE atom, ONLY: nlcc
|
2003-11-09 18:42:50 +08:00
|
|
|
USE wavefunctions_module, ONLY: evc
|
2003-12-11 19:10:03 +08:00
|
|
|
use io_files, only: prefix, iunpun, iunres
|
2003-01-20 05:58:50 +08:00
|
|
|
use cgcom
|
|
|
|
use funct
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
!
|
|
|
|
integer :: i, l, nt, kpoint
|
|
|
|
logical :: exst
|
2003-04-13 03:25:08 +08:00
|
|
|
character (len=256) :: filint
|
2003-01-20 05:58:50 +08:00
|
|
|
real(kind=DP) :: rhotot, dmxc
|
|
|
|
external dmxc
|
|
|
|
!
|
|
|
|
call start_clock('cg_setup')
|
|
|
|
!
|
|
|
|
! convert masses to atomic units
|
|
|
|
!
|
|
|
|
call DSCAL(ntyp,amconv,amass,1)
|
|
|
|
!
|
|
|
|
! sum self-consistent part (vr) and local part (vltot) of potential
|
|
|
|
!
|
|
|
|
call set_vrs(vrs,vltot,vr,nrxx,nspin,doublegrid)
|
|
|
|
!
|
|
|
|
! allocate memory for various arrays
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate (dmuxc( nrxx))
|
|
|
|
allocate (dvpsi( npwx, nbnd))
|
|
|
|
allocate ( dpsi( npwx, nbnd))
|
|
|
|
allocate ( auxr( nrxx))
|
|
|
|
allocate ( aux2( nrxx))
|
|
|
|
allocate ( aux3( nrxx))
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! allocate memory for gradient corrections (if needed)
|
|
|
|
!
|
|
|
|
if (igcx.ne.0 .or. igcc.ne.0) then
|
2003-02-08 00:04:36 +08:00
|
|
|
allocate ( dvxc_rr(nrxx,nspin,nspin))
|
|
|
|
allocate ( dvxc_sr(nrxx,nspin,nspin))
|
|
|
|
allocate ( dvxc_ss(nrxx,nspin,nspin))
|
|
|
|
allocate ( dvxc_s (nrxx,nspin,nspin))
|
|
|
|
allocate ( grho (3, nrxx, nspin))
|
2003-01-20 05:58:50 +08:00
|
|
|
end if
|
|
|
|
!
|
|
|
|
!
|
|
|
|
! initialize structure factor array
|
|
|
|
!
|
|
|
|
call struc_fact ( nat, tau, ntyp, ityp, ngm, g, bg, &
|
|
|
|
& nr1, nr2, nr3, strf, eigts1, eigts2, eigts3 )
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
! compute drhocore/dtau for each atom type (if needed)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
nlcc_any = .false.
|
|
|
|
do nt=1,ntyp
|
|
|
|
nlcc_any = nlcc_any .or. nlcc(nt)
|
|
|
|
end do
|
|
|
|
!!! if (nlcc_any) call set_drhoc(xq)
|
|
|
|
!
|
|
|
|
! local potential
|
|
|
|
!
|
|
|
|
call init_vloc
|
|
|
|
!
|
|
|
|
call init_us_1
|
|
|
|
!
|
|
|
|
call newd
|
|
|
|
!
|
|
|
|
! derivative of the xc potential
|
|
|
|
!
|
2003-04-13 03:25:08 +08:00
|
|
|
dmuxc(:) = 0.d0
|
2003-01-20 05:58:50 +08:00
|
|
|
do i = 1,nrxx
|
|
|
|
rhotot = rho(i,current_spin)+rho_core(i)
|
|
|
|
if ( rhotot.gt. 1.d-30 ) dmuxc(i)= dmxc( rhotot)
|
|
|
|
if ( rhotot.lt.-1.d-30 ) dmuxc(i)=-dmxc(-rhotot)
|
|
|
|
end do
|
|
|
|
!
|
|
|
|
! initialize data needed for gradient corrections
|
|
|
|
!
|
|
|
|
call cg_setupdgc
|
|
|
|
!
|
|
|
|
iunres=88
|
|
|
|
!
|
|
|
|
! open the wavefunction file (already existing)
|
|
|
|
!
|
|
|
|
lrwfc=2*nbnd*npwx
|
2003-02-08 00:04:36 +08:00
|
|
|
filint = trim(prefix) //'.wfc'
|
2003-01-20 05:58:50 +08:00
|
|
|
call diropn(iunpun,filint,lrwfc,exst)
|
2003-12-11 19:10:03 +08:00
|
|
|
if(.not.exst) call errore('main','file '//filint//' not found',1)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! read wave functions and calculate indices
|
|
|
|
!
|
|
|
|
kpoint=1
|
|
|
|
call davcio(evc,lrwfc,iunpun,kpoint,-1)
|
|
|
|
close(unit=iunpun,status='keep')
|
|
|
|
call gk_sort (xk(1,kpoint),ngm,g,ecutwfc/tpiba2,npw,igk,g2kin)
|
|
|
|
!
|
|
|
|
! Kleinman-Bylander PPs
|
|
|
|
!
|
|
|
|
call init_us_2 (npw, igk, xk(1,kpoint), vkb)
|
|
|
|
!
|
|
|
|
call stop_clock('cg_setup')
|
|
|
|
!
|
|
|
|
return
|
|
|
|
end subroutine cg_setup
|