! ! Copyright (C) 2002-2008 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 pseudo_types ! this module contains the definitions of several TYPE structures, ! together with their allocation/deallocation routines USE kinds, ONLY: DP use radial_grids, ONLY: radial_grid_type IMPLICIT NONE SAVE ! ! Additional data to make a PAW setup out of an US pseudo, ! they are all stored on a radial grid: TYPE paw_in_upf REAL(DP),POINTER :: ae_rho_atc(:) ! AE core charge (pseudo ccharge ! is already included in upf) REAL(DP),POINTER :: pfunc(:,:,:),&! Psi_i(r)*Psi_j(r) ptfunc(:,:,:) ! as above, but for pseudo REAL(DP),POINTER :: ae_vloc(:) ! AE local potential (pseudo vloc ! is already included in upf) REAL(DP),POINTER :: oc(:) ! starting occupation used to init becsum ! they differ from US ones because they ! are indexed on BETA functions, non on WFC REAL(DP),POINTER :: augmom(:,:,:) ! multipole AE-pseudo (i,j,l=0:2*lmax) REAL(DP) :: raug ! augfunction max radius INTEGER :: iraug ! index on rgrid closer to, and >, raug INTEGER :: lmax_aug ! max angmom of augmentation functions, it is == ! to 2* max{l of pseudized wavefunctions} ! note that nqlc of upf also includes the angmom of ! empty virtual channel used to generate local potential REAL(DP) :: core_energy ! constant to add in order to get all-electron energy CHARACTER(len=12):: augshape ! shape of augmentation charge END TYPE paw_in_upf TYPE pseudo_upf CHARACTER(LEN=80):: generated='' ! generator software CHARACTER(LEN=80):: author='' ! pseudopotential's author CHARACTER(LEN=80):: date='' ! generation date CHARACTER(LEN=80):: comment='' ! author's comment CHARACTER(LEN=2) :: psd='' ! Element label CHARACTER(LEN=20) :: typ='' ! Pseudo type ( NC or US or PAW) CHARACTER(len=6) :: rel='' ! relativistic: {no|scalar|full} LOGICAL :: tvanp ! .true. if Ultrasoft LOGICAL :: tcoulombp ! .true. if Coulomb 1/r potential LOGICAL :: nlcc ! Non linear core corrections CHARACTER(LEN=20) :: dft ! Exch-Corr type REAL(DP) :: zp ! z valence REAL(DP) :: etotps ! total energy REAL(DP) :: ecutwfc ! suggested cut-off for wfc REAL(DP) :: ecutrho ! suggested cut-off for rho ! CHARACTER(len=11) :: nv ! UPF file three-digit version i.e. 2.0.0 INTEGER :: lmax ! maximum l component in beta INTEGER :: lmax_rho ! max l componet in charge (should be 2*lmax) ! Wavefunctions and projectors INTEGER :: nwfc ! number of atomic wavefunctions INTEGER :: nbeta ! number of projectors INTEGER, POINTER :: kbeta(:) ! kbeta(nbeta) see below INTEGER :: kkbeta ! kkbeta=max(kbeta(:)) ! kbeta<=mesh is the number of grid points for each beta function ! beta(r,nb) = 0 for r > r(kbeta(nb)) ! kkbeta<=mesh is the largest of such number so that for all beta ! beta(r,nb) = 0 for r > r(kkbeta) ! INTEGER, POINTER :: lll(:) ! lll(nbeta) l of each projector REAL(DP), POINTER :: beta(:,:) ! beta(mesh,nbeta) projectors ! CHARACTER(LEN=2), POINTER :: els(:) ! els(nwfc) label of wfc CHARACTER(LEN=2), POINTER :: els_beta(:) ! els(nbeta) label of beta INTEGER, POINTER :: nchi(:) ! lchi(nwfc) value of pseudo-n for wavefcts INTEGER, POINTER :: lchi(:) ! lchi(nwfc) value of l for wavefcts REAL(DP), POINTER :: oc(:) ! oc(nwfc) occupancies for wavefcts REAL(DP), POINTER :: epseu(:) ! pseudo one-particle energy (nwfc) REAL(DP), POINTER :: rcut_chi(:)! rcut_chi(nwfc) cutoff innner radius REAL(DP), POINTER :: rcutus_chi(:)! rcutus_chi(nwfc) ultrasoft outer radius ! Chi and rho_at are only used for initial density and initial wfcs: REAL(DP), POINTER :: chi(:,:) ! chi(mesh,nwfc) atomic wavefcts REAL(DP), POINTER :: rho_at(:) ! rho_at(mesh) atomic charge ! Minimal radial grid: INTEGER :: mesh ! number of points in the radial mesh REAL(DP) :: xmin ! the minimum x of the linear mesh REAL(DP) :: rmax ! the maximum radius of the mesh REAL(DP) :: zmesh ! the nuclear charge used for mesh REAL(DP) :: dx ! the deltax of the linear mesh REAL(DP), POINTER :: r(:) ! r(mesh) radial grid REAL(DP), POINTER :: rab(:) ! rab(mesh) dr(x)/dx (x=linear grid) ! Pseudized core charge REAL(DP), POINTER :: rho_atc(:) ! rho_atc(mesh) atomic core charge ! Local potential INTEGER :: lloc ! L of channel used to generate local potential ! (if < 0 it was generated by smoothing AE potential) REAL(DP) :: rcloc ! vloc = v_ae for r > rcloc REAL(DP), POINTER :: vloc(:) ! vloc(mesh) local atomic potential ! REAL(DP), POINTER :: dion(:,:) ! dion(nbeta,nbeta) atomic D_{mu,nu} ! Augmentation LOGICAL :: q_with_l ! if .true. qfunc is pseudized in ! different ways for different l INTEGER :: nqf ! number of Q coefficients INTEGER :: nqlc ! number of angular momenta in Q REAL(DP):: qqq_eps ! qfunc is null if its norm is .lt. qqq_eps REAL(DP), POINTER :: rinner(:) ! rinner(0:2*lmax) r_L REAL(DP), POINTER :: qqq(:,:) ! qqq(nbeta,nbeta) q_{mu,nu} ! Augmentation without L dependecy REAL(DP), POINTER :: qfunc(:,:) ! qfunc(mesh,nbeta*(nbeta+1)/2) ! Q_{mu,nu}(|r|) function for |r|> r_L ! Augmentation depending on L (optional, compulsory for PAW) REAL(DP), POINTER :: qfuncl(:,:,:)! qfuncl(mesh,nbeta*(nbeta+1)/2,l) ! Q_{mu,nu}(|r|) function for |r|> r_L ! Analitycal coeffs cor small r expansion of qfunc (Vanderbilt's code) REAL(DP), POINTER :: qfcoef(:,:,:,:) ! qfcoef(nqf,0:2*lmax,nbeta,nbeta) ! coefficients for Q for |r|