quantum-espresso/Modules/uspp.f90

342 lines
11 KiB
Fortran
Raw Normal View History

!
! Copyright (C) 2004-2011 Quantum ESPRESSO 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 .
!
MODULE uspp_param
!
! ... Ultrasoft and Norm-Conserving pseudopotential parameters
!
SCF with Projector-Augmented Wave Pseudopotential (PAW) routines added. This means that a lot of routines have been modified and a few files have been added. During the year several people have contributed to this code, mainly Guido Fratesi, Ricardo Mazzarello, Stefano de Gironcoli, Andrea Dal Corso and me (Lorenzo Paulatto). A brief report of modified or added files follows, further down you will find a loger report of modifications that was necessary to merge develop_PAW branch with the current CVS version. Current version is not 100% functional, but it doesn't brake anything else and can be used to generate and test PAW pseudopotential. ************************************* *** Brief report of modifications *** ************************************* Modified files: PW/clean_pw.f90 PW/electrons.f90 PW/print_clock_pw.f90 PW/hinit0.f90 PW/potinit.f90 PW/newd.f90 PW/summary.f90 PW/setup.f90 PW/read_pseudo.f90 PW/init_us_1.f90 PW/init_run.f90 PW/mix_rho.f90 atomic/atomic_paw.f90 atomic/write_paw_recon.f90 atomic/ld1_writeout.f90 atomic/write_resultsps.f90 atomic/ld1inc.f90 atomic/ld1_readin.f90 atomic/gener_pseudo.f90 atomic/parameters.f900 atomic/run_pseudo.f900 atomic/set_rho_core.f90 atomic/pseudovloc.f90 Modules/read_upf.f90 Modules/uspp.f90 Modules/pseudo_types.f90 Modules/parameters.f90 Added files: PW/grid_paw_routines.f90 PW/rad_paw_routines.f90 Modules/grid_paw_variables.f90 Modules/read_paw.f90 Added files that will be removed: PW/rad_paw_trash.f90 PW/paw_xc.f90 Examples: examples/PAWexample contains a full test of PAW pseudopotential for Oxygen. The test consist in these tasks: - 2 norm conserving, 2 US and 4 PAW pseudopotentials are generated and tested in ld1 - pw test for an isolated O atom at different cutoffs - pw test for an O2 molecule at different O-O distance please read examples/PAWexample/README for (a few) details. NOTES: 1. new modifications to atomic_paw (and related) from ADC have been rolled back, as they were breaking a lot of things, I will reintroduce them later when I am sure that everything works properly. 2. the files PW/paw_xc.f90 and Modules/rad_paw_trash.f90 will be removed in the next few weeks. TODO: 1. use new ld1 XC code as much as possible, and remove legacy XC routines from rad_paw_routines 2. full self-consistency with radial energies 3. make new Harris-Foulkes estimate paw-aware 4. provide some kind error estimate 5. FORCES and stress!! (require symmetrization of becsums) 6. cleanup ************************ *** merge of PW code *** ************************ Versions notation: OLD=version from 2 years ago used as reference to generate the patches NEW=CURRENT=current trunk version PAW=current develop_PAW version Note: pseudo-potential input and allocation routines changed a lot in the last years, this is a diagram: OLD:PW/readin ~~> PAW:PW/read_pseudo --> disappears pops out --> PAW:PW/readin ~~> NEW:PW/read_pseudo added files: Modules/read_paw.f90 (contains module read_paw_module with subroutines paw_io nullify_pseudo_paw, allocate_pseudo_paw and deallocate_pseudo_paw previously in removed file Modules/readpseudo.f90. Also contains module paw_to_internal with subroutine set_pseudo_paw, previously in upf_to_internal.f90) PW/paw_xc.f90 (contains OLD=PAW xc and gcxc routines as adapting paw grid code to use new routines was very error prone and quite worthless, as it has to be removed anyway) Conflicts reported by CVS during merge: DONE */Makefiles (all replaced with new, redone by hand) DONE flib/functionals.f90 (nothing to do) DONE Modules/functionals.f90 (RNV == replaced with NEW version) DONE Modules/atom.f90 (trivial: duped rgrid) DONE Modules/autopilot.f90 (trivial) DONE Modules/bfgs_module.f90 (RNV) DONE Modules/cell_base.f90 (RNV) DONE Modules/check_stop.f90 (RNV) DONE Modules/constants.f90 (RNV) DONE Modules/constraints_module.f90 (RNV) DONE Modules/energies.f90 (RNV) DONE Modules/input_parameters.f90 (RNV) DONE Modules/ions_base.f90 (RNV, has 3 new subs) DONE Modules/ions_nose.f90 (RNV) DONE Modules/parameters.f90 (actually RNV) DONE Modules/path_base.f90 (RNV) DONE Modules/path_opt_routines.f90 (RNV) DONE Modules/path_reparametrisation.f90 (RNV) DONE Modules/path_variables.f90 (RNV) DONE Modules/pseudo_types.f90 (cleaned double def of paw_t) DONE Modules/read_cards.f90 (RNV) DONE Modules/read_namelists.f90 (checked and RNV) DONE Modules/uspp.f90 (trivial) DONE Modules/xml_io_base.f90 (RNV) DONE PW/read_pseudo.f90 (merged by hand with PAW PW/readin) DONE PW/bp_calc_btq.f90 (trivial) DONE PW/c_bands.f90 (actually RNV) DONE PW/ccgdiagg.f90 (RNV) DONE PW/cegterg.f90 (RNV) DONE PW/cft3s.f90 (RNV) DONE PW/cinitcgg.f90 (RNV) DONE PW/c_phase_field.f90 (RNV) DONE PW/divide_et_impera.f90 (nothing to do?) DONE PW/exx.f90 (RNV) DONE PW/hinit0.f90 (easy) DONE PW/h_psi.f90 (RNV) DONE PW/init_run.f90 (easy) DONE PW/kpoint_grid.f90 (nothing to do?) DONE PW/newd.f90 (required mod in newd_paw_grid, CHECK!!) DONE PW/openfil.f90 (actually RNV) DONE PW/paw.f90 (actually RNV) DONE PW/punch.f90 (RNV) DONE PW/pwscf.f90 (quite RNV) DONE PW/set_kup_and_kdw.f90 (RNV) DONE PW/setup.f90 (RNV + 2 line merged by hand) DONE PW/sgama.f9 (actually RNV) DONE PW/sgam_at_mag.f90 (actually RNV) DONE PW/stop_run.f90 (actually RNV) DONE PW/stres_gradcorr.f90 (actually RNV) DONE PW/symrho_mag.f90 (nothing to do) DONE PW/v_of_rho.f90 (RNV) DONE PW/compute_fes_grads.f90 (RNV) DONE PW/gradcorr.f90 (RNV) DONE PW/input.f90 (RNV) DONE PW/pw_restart.f90 (RNV) DONE PW/read_ncpp.f90 (actually RNV) DONE PW/summary.f90 (RNV + inserted new PP type) DONE PW/wfcinit.f90 (RNV) the hard ones: DONE PW/electrons.f90 (adapted code to new syntaxes, a lot of cleanup, removed some PAW junk that can be readded later, removed parts that were applyed twice, or had been removed in trunk, the rhog allocations and usage may need fixes) DONE PW/mix_rho.f90 (merged tauk and paw additions, a bit of cleanup and smarter variables names) DONE PW/init_us_1.f90 (qtot redefined with "triangular" index nb,mb-->ijv) modified for compiling: Modules/io_files.f90 (depatched) PW/pwcom.f90 (depatched) Modules/parameters.f90 (temporary readded cp_lmax = lmaxx+1) PW/newd.f90 (merge was wrong, redone mostly by hand) PW/read_ncpp.f90 (depatched) PW/read_pseudo (small fixes) PW/sgam_at_mag.f90 (depatched) PW/sgama.f90 (depatched) PW/stres_gradcorr.f90 (depatched) modified for running: PW/clean_pw.f90 (added call to deallocate_paw_internals) Modifications to PAW routines: 1. compute_onecenter_charges and compute_onecenter_charges modified to comply with new structure of v_xc (in v_of_rho.f90), requiring new g-space densities to be saved and computed --> using old xc routines as this code will be removed. 2. qrad size has changed, prad and ptrad had to be changed accordingly. 3. several minor modifications to use new radial grid structure. 4. infomsg arguments changed, very funny bug followed. 5. added new routine deallocate_paw_internals, called by PW/clean_pw.f90 required to run pp.x with more than one q-point(and good programming practice) ************************* *** merge of LD1 code *** ************************* 2nd try: atomic code replaced with current version, then merge by hand the files that are used by paw subsystem: * atomic_paw.f90 (replaced with most recent version from develop_PAW routine us2paw and paw2us taken from newer trunk version, a lot of minor changes.) * gener_pseudo.f90 (fixes) * ld1inc.f90 (PAW variables added) * ld1_readin.f90 (PAW variables added, I am not sure if lpaw should go in input or inputp namelist) * ld1_writeout.f90 (it was only necessary to add a few lines) * pseudovloc.f90 (nothing to do) * run_pseudo.f90 (almost nothing to do) * set_rho_core.f90 (readded a few lines for lnc2paw) * write_paw_recon.f90 (nothing to do) * write_resultsps.f90 (nothing to do: trunk version is more PAW-aware than PAW version) Main problems were found in subroutines run_pseudo and gen_pseudo, a little code had to be rewritten to comply with new variable names and fix with merge. TODO: fix atomic_paw routines to use minimal allocated arrays insetad of ndmx sized ones; try to use the pawet as much as possible. Remove test lines and other garbage. Find a fix for PAW2. The first week of september Andrea Dal Corso uploaded a few modifications to the atomic_paw routines. I had to rollback them as the structure of atomic_paw has changed a lot and reimplementing them is probably easier and definitely safer than fixing everything. I will do it soon, I swear! LP git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@4257 c92efa57-630b-4861-b058-cf58834340f0
2007-09-18 18:05:46 +08:00
USE kinds, ONLY : DP
USE parameters, ONLY : npsx
USE pseudo_types, ONLY : pseudo_upf
!
SAVE
PUBLIC :: n_atom_wfc
!
TYPE (pseudo_upf), ALLOCATABLE, TARGET :: upf(:)
INTEGER :: &
nh(npsx), &! number of beta functions per atomic type
nhm, &! max number of different beta functions per atom
nbetam, &! max number of beta functions
iver(3,npsx) ! version of the atomic code
INTEGER :: &
lmaxkb, &! max angular momentum
lmaxq ! max angular momentum + 1 for Q functions
LOGICAL :: &
newpseudo(npsx), &! if .TRUE. multiple projectors are allowed
oldvan(npsx) ! old version of Vanderbilt PPs, using
! Herman-Skillman grid - obsolescent
INTEGER :: &
nvb, &! number of species with Vanderbilt PPs (CPV)
ish(npsx) ! for each specie the index of the first beta
! function: ish(1)=1, ish(i)=1+SUM(nh(1:i-1))
CONTAINS
!
!----------------------------------------------------------------------------
FUNCTION n_atom_wfc( nat, ityp, noncolin )
!----------------------------------------------------------------------------
!
! ... Find number of starting atomic orbitals
!
IMPLICIT NONE
!
INTEGER, INTENT(IN) :: nat, ityp(nat)
LOGICAL, INTENT(IN), OPTIONAL :: noncolin
INTEGER :: n_atom_wfc
!
INTEGER :: na, nt, n
LOGICAL :: non_col
!
!
non_col = .FALSE.
IF ( PRESENT (noncolin) ) non_col=noncolin
n_atom_wfc = 0
!
DO na = 1, nat
!
nt = ityp(na)
!
DO n = 1, upf(nt)%nwfc
!
IF ( upf(nt)%oc(n) >= 0.D0 ) THEN
!
IF ( non_col ) THEN
!
IF ( upf(nt)%has_so ) THEN
!
n_atom_wfc = n_atom_wfc + 2 * upf(nt)%lchi(n)
!
IF ( ABS( upf(nt)%jchi(n)-upf(nt)%lchi(n) - 0.5D0 ) < 1.D-6 ) &
n_atom_wfc = n_atom_wfc + 2
!
ELSE
!
n_atom_wfc = n_atom_wfc + 2 * ( 2 * upf(nt)%lchi(n) + 1 )
!
END IF
!
ELSE
!
n_atom_wfc = n_atom_wfc + 2 * upf(nt)%lchi(n) + 1
!
END IF
END IF
END DO
END DO
!
RETURN
!
END FUNCTION n_atom_wfc
END MODULE uspp_param
! <<<<<<<<<<<<<<<~~~~<<<<<<<<<<<<<<<<-----------------
MODULE uspp
!
! Ultrasoft PPs:
! - Clebsch-Gordan coefficients "ap", auxiliary variables "lpx", "lpl"
! - beta and q functions of the solid
!
USE kinds, ONLY: DP
USE parameters, ONLY: lmaxx, lqmax
IMPLICIT NONE
PRIVATE
SAVE
PUBLIC :: nlx, lpx, lpl, ap, aainit, indv, nhtol, nhtolm, indv_ijkb0, &
nkb, nkbus, vkb, dvan, deeq, qq, nhtoj, ijtoh, beta, &
becsum, deallocate_uspp
PUBLIC :: okvan, nlcc_any
PUBLIC :: qq_so, dvan_so, deeq_nc
PUBLIC :: dbeta
INTEGER, PARAMETER :: &
nlx = (lmaxx+1)**2, &! maximum number of combined angular momentum
mx = 2*lqmax-1 ! maximum magnetic angular momentum of Q
INTEGER :: &! for each pair of combined momenta lm(1),lm(2):
lpx(nlx,nlx), &! maximum combined angular momentum LM
lpl(nlx,nlx,mx) ! list of combined angular momenta LM
REAL(DP) :: &
ap(lqmax*lqmax,nlx,nlx)
! Clebsch-Gordan coefficients for spherical harmonics
!
INTEGER :: nkb, &! total number of beta functions, with struct.fact.
nkbus ! as above, for US-PP only
!
INTEGER, ALLOCATABLE ::&
indv(:,:), &! indes linking atomic beta's to beta's in the solid
nhtol(:,:), &! correspondence n <-> angular momentum l
nhtolm(:,:), &! correspondence n <-> combined lm index for (l,m)
ijtoh(:,:,:), &! correspondence beta indexes ih,jh -> composite index ijh
indv_ijkb0(:) ! first beta (index in the solid) for each atom
!
LOGICAL :: &
okvan = .FALSE.,& ! if .TRUE. at least one pseudo is Vanderbilt
nlcc_any=.FALSE. ! if .TRUE. at least one pseudo has core corrections
!
COMPLEX(DP), ALLOCATABLE, TARGET :: &
vkb(:,:) ! all beta functions in reciprocal space
REAL(DP), ALLOCATABLE :: &
becsum(:,:,:) ! \sum_i f(i) <psi(i)|beta_l><beta_m|psi(i)>
REAL(DP), ALLOCATABLE :: &
dvan(:,:,:), &! the D functions of the solid
deeq(:,:,:,:), &! the integral of V_eff and Q_{nm}
qq(:,:,:), &! the q functions in the solid
nhtoj(:,:) ! correspondence n <-> total angular momentum
!
COMPLEX(DP), ALLOCATABLE :: & ! variables for spin-orbit/noncolinear case:
qq_so(:,:,:,:), &! Q_{nm}
dvan_so(:,:,:,:), &! D_{nm}
deeq_nc(:,:,:,:) ! \int V_{eff}(r) Q_{nm}(r) dr
!
! spin-orbit coupling: qq and dvan are complex, qq has additional spin index
! noncolinear magnetism: deeq is complex (even in absence of spin-orbit)
!
REAL(DP), ALLOCATABLE :: &
beta(:,:,:) ! beta functions for CP (without struct.factor)
REAL(DP), ALLOCATABLE :: &
dbeta(:,:,:,:,:) ! derivative of beta functions w.r.t. cell for CP (without struct.factor)
!
CONTAINS
!
!-----------------------------------------------------------------------
subroutine aainit(lli)
!-----------------------------------------------------------------------
!
! this routine computes the coefficients of the expansion of the product
! of two real spherical harmonics into real spherical harmonics.
!
! Y_limi(r) * Y_ljmj(r) = \sum_LM ap(LM,limi,ljmj) Y_LM(r)
!
! On output:
! ap the expansion coefficients
! lpx for each input limi,ljmj is the number of LM in the sum
! lpl for each input limi,ljmj points to the allowed LM
!
! The indices limi,ljmj and LM assume the order for real spherical
! harmonics given in routine ylmr2
!
USE matrix_inversion
implicit none
!
! input: the maximum li considered
!
integer :: lli
!
! local variables
!
integer :: llx, l, li, lj
real(DP) , allocatable :: r(:,:), rr(:), ylm(:,:), mly(:,:)
! an array of random vectors: r(3,llx)
! the norm of r: rr(llx)
! the real spherical harmonics for array r: ylm(llx,llx)
! the inverse of ylm considered as a matrix: mly(llx,llx)
!
if (lli < 0) call errore('aainit','lli not allowed',lli)
if (lli*lli > nlx) call errore('aainit','nlx is too small ',lli*lli)
llx = (2*lli-1)**2
if (2*lli-1 > lqmax) &
call errore('aainit','ap leading dimension is too small',llx)
allocate (r( 3, llx ))
allocate (rr( llx ))
allocate (ylm( llx, llx ))
allocate (mly( llx, llx ))
r(:,:) = 0.0_DP
ylm(:,:) = 0.0_DP
mly(:,:) = 0.0_DP
ap(:,:,:)= 0.0_DP
! - generate an array of random vectors (uniform deviate on unitary sphere)
call gen_rndm_r(llx,r,rr)
! - generate the real spherical harmonics for the array: ylm(ir,lm)
call ylmr2(llx,llx,r,rr,ylm)
!- store the inverse of ylm(ir,lm) in mly(lm,ir)
call invmat(llx, ylm, mly)
!- for each li,lj compute ap(l,li,lj) and the indices, lpx and lpl
do li = 1, lli*lli
do lj = 1, lli*lli
lpx(li,lj)=0
do l = 1, llx
ap(l,li,lj) = compute_ap(l,li,lj,llx,ylm,mly)
if (abs(ap(l,li,lj)) > 1.d-3) then
lpx(li,lj) = lpx(li,lj) + 1
if (lpx(li,lj) > mx) &
call errore('aainit','mx dimension too small', lpx(li,lj))
lpl(li,lj,lpx(li,lj)) = l
end if
end do
end do
end do
deallocate(mly)
deallocate(ylm)
deallocate(rr)
deallocate(r)
return
end subroutine aainit
!
!-----------------------------------------------------------------------
subroutine gen_rndm_r(llx,r,rr)
!-----------------------------------------------------------------------
! - generate an array of random vectors (uniform deviate on unitary sphere)
!
USE constants, ONLY: tpi
USE random_numbers, ONLY: randy
implicit none
!
! first the I/O variables
!
integer :: llx ! input: the dimension of r and rr
real(DP) :: &
r(3,llx), &! output: an array of random vectors
rr(llx) ! output: the norm of r
!
! here the local variables
!
integer :: ir
real(DP) :: costheta, sintheta, phi
do ir = 1, llx
costheta = 2.0_DP * randy() - 1.0_DP
sintheta = SQRT ( 1.0_DP - costheta*costheta)
phi = tpi * randy()
r (1,ir) = sintheta * cos(phi)
r (2,ir) = sintheta * sin(phi)
r (3,ir) = costheta
rr(ir) = 1.0_DP
end do
return
end subroutine gen_rndm_r
!
!-----------------------------------------------------------------------
function compute_ap(l,li,lj,llx,ylm,mly)
!-----------------------------------------------------------------------
!- given an l and a li,lj pair compute ap(l,li,lj)
implicit none
!
! first the I/O variables
!
integer :: &
llx, &! the dimension of ylm and mly
l,li,lj ! the arguments of the array ap
real(DP) :: &
compute_ap, &! this function
ylm(llx,llx),&! the real spherical harmonics for array r
mly(llx,llx) ! the inverse of ylm considered as a matrix
!
! here the local variables
!
integer :: ir
compute_ap = 0.0_DP
do ir = 1,llx
compute_ap = compute_ap + mly(l,ir)*ylm(ir,li)*ylm(ir,lj)
end do
return
end function compute_ap
!
!-----------------------------------------------------------------------
SUBROUTINE deallocate_uspp()
!-----------------------------------------------------------------------
!
IF( ALLOCATED( nhtol ) ) DEALLOCATE( nhtol )
IF( ALLOCATED( indv ) ) DEALLOCATE( indv )
IF( ALLOCATED( nhtolm ) ) DEALLOCATE( nhtolm )
IF( ALLOCATED( nhtoj ) ) DEALLOCATE( nhtoj )
IF( ALLOCATED( indv_ijkb0 ) ) DEALLOCATE( indv_ijkb0 )
IF( ALLOCATED( ijtoh ) ) DEALLOCATE( ijtoh )
IF( ALLOCATED( vkb ) ) DEALLOCATE( vkb )
IF( ALLOCATED( becsum ) ) DEALLOCATE( becsum )
IF( ALLOCATED( qq ) ) DEALLOCATE( qq )
IF( ALLOCATED( dvan ) ) DEALLOCATE( dvan )
IF( ALLOCATED( deeq ) ) DEALLOCATE( deeq )
IF( ALLOCATED( qq_so ) ) DEALLOCATE( qq_so )
IF( ALLOCATED( dvan_so ) ) DEALLOCATE( dvan_so )
IF( ALLOCATED( deeq_nc ) ) DEALLOCATE( deeq_nc )
IF( ALLOCATED( beta ) ) DEALLOCATE( beta )
IF( ALLOCATED( dbeta ) ) DEALLOCATE( dbeta )
!
END SUBROUTINE deallocate_uspp
!
END MODULE uspp