! ! Copyright (C) 2002 CP90 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 bhs ! analytical BHS pseudopotential parameters use parameters, only: nsx implicit none save real(kind=8) :: rc1(nsx), rc2(nsx), wrc1(nsx), wrc2(nsx), & rcl(3,nsx,3), al(3,nsx,3), bl(3,nsx,3) integer :: lloc(nsx) end module bhs module core implicit none save ! nlcc_any = 0 no core correction on any atom ! rhocb = core charge in G space (box grid) logical :: nlcc_any real(kind=8), allocatable:: rhocb(:,:) contains subroutine deallocate_core() IF( ALLOCATED( rhocb ) ) DEALLOCATE( rhocb ) end subroutine end module core module cvan ! ionic pseudo-potential variables use parameters, only: nsx implicit none save ! nvb = number of species with Vanderbilt PPs ! nh(is) = number of beta functions, including Y_lm, for species is ! ish(is)= used for indexing the nonlocal projectors betae ! with contiguous indices inl=ish(is)+(iv-1)*na(is)+1 ! where "is" is the species and iv=1,nh(is) ! nhx = max value of nh(np) ! nhsavb = total number of Vanderbilt nonlocal projectors ! nhsa = total number of nonlocal projectors for all atoms integer nvb, nhsavb, ish(nsx), nh(nsx), nhsa, nhx ! indv : indv( ind,is)=beta function (without Y_lm) for projector ind ! nhtol : nhtol( ind,is)=value of l for projector ind of species is ! nhtolm: nhtolm(ind,is)=cobined lm index in Y_lm for projector ind integer, allocatable:: nhtol(:,:), indv(:,:), nhtolm(:,:) ! beta = nonlocal projectors in g space without e^(-ig.r) factor ! qq = ionic Q_ij for each species (Vanderbilt only) ! dvan = ionic D_ij for each species (Vanderbilt only) real(kind=8), allocatable:: beta(:,:,:), qq(:,:,:), dvan(:,:,:) contains subroutine deallocate_cvan() IF( ALLOCATED( nhtol ) ) DEALLOCATE( nhtol ) IF( ALLOCATED( indv ) ) DEALLOCATE( indv ) IF( ALLOCATED( nhtolm ) ) DEALLOCATE( nhtolm ) IF( ALLOCATED( beta ) ) DEALLOCATE( beta ) IF( ALLOCATED( qq ) ) DEALLOCATE( qq ) IF( ALLOCATED( dvan ) ) DEALLOCATE( dvan ) end subroutine end module cvan module elct use electrons_base, only: nspin, nel, nupdwn, iupdwn use electrons_base, only: n => nbnd, nx => nbndx implicit none save ! f = occupation numbers ! qbac = background neutralizing charge real(kind=8), allocatable:: f(:) real(kind=8) qbac ! nspin = number of spins (1=no spin, 2=LSDA) ! nel(nspin) = number of electrons (up, down) ! nupdwn= number of states with spin up (1) and down (2) ! iupdwn= first state with spin (1) and down (2) ! n = total number of electronic states ! nx = if n is even, nx=n ; if it is odd, nx=n+1 ! nx is used only to dimension arrays ! ispin = spin of each state integer, allocatable:: ispin(:) ! contains subroutine deallocate_elct() IF( ALLOCATED( f ) ) DEALLOCATE( f ) IF( ALLOCATED( ispin ) ) DEALLOCATE( ispin ) return end subroutine ! end module elct module gvec use cell_base, only: tpiba, tpiba2 use reciprocal_vectors, only: & gl, g, gx, g2_g, mill_g, mill_l, ig_l2g, igl, bi1, bi2, bi3 use reciprocal_vectors, only: deallocate_recvecs use recvecs_indexes, only: np, nm, in1p, in2p, in3p, deallocate_recvecs_indexes use gvecp, only: & ng => ngm, & ngl => ngml, & ng_g => ngmt ! tpiba = 2*pi/alat ! tpiba2 = (2*pi/alat)**2 ! ng = number of G vectors for density and potential ! ngl = number of shells of G ! G-vector quantities for the thick grid - see also doc in ggen ! g = G^2 in increasing order (in units of tpiba2=(2pi/a)^2) ! gl = shells of G^2 ( " " " " " ) ! gx = G-vectors ( " " " tpiba =(2pi/a) ) ! ! g2_g = all G^2 in increasing order, replicated on all procs ! mill_g = miller index of G vecs (increasing order), replicated on all procs ! mill_l = miller index of G vecs local to the processors ! ig_l2g = "l2g" means local to global, this array convert a local ! G-vector index into the global index, in other words ! the index of the G-v. in the overall array of G-vectors ! bi? = base vector used to generate the reciprocal space ! ! np = fft index for G> ! nm = fft index for G< ! in1p,in2p,in3p = G components in crystal axis ! implicit none save contains subroutine deallocate_gvec CALL deallocate_recvecs( ) CALL deallocate_recvecs_indexes( ) end subroutine end module gvec module ncprm use parameters, only: nsx, ndmx, nqfx, nbrx, lqmax implicit none save ! ! lqmax: maximum angular momentum of Q (Vanderbilt augmentation charges) ! nqfx : maximum number of coefficients in Q smoothing ! nbrx : maximum number of distinct radial beta functions ! ndmx: maximum number of points in the radial grid ! ! nbeta number of beta functions (sum over all l) ! kkbeta last radial mesh point used to describe functions ! which vanish outside core ! nqf coefficients in Q smoothing ! nqlc angular momenta present in Q smoothing ! lll lll(j) is l quantum number of j'th beta function ! lqx highest angular momentum that is present in Q functions ! lmaxkb highest angular momentum that is present in beta functions integer :: nbeta(nsx), kkbeta(nsx), & nqf(nsx), nqlc(nsx), lll(nbrx,nsx), lqx, lmaxkb ! dion bare pseudopotential D_{\mu,\nu} parameters ! (ionic and screening parts subtracted out) ! betar the beta function on a r grid (actually, r*beta) ! qqq Q_ij matrix ! qfunc Q_ij(r) function (for r>rinner) ! rinner radius at which to cut off partial core or Q_ij ! ! qfcoef coefficients to pseudize qfunc for different total ! angular momentum (for r