TABLE OF CONTENTS


ABINIT/defs_datatypes [ Modules ]

[ Top ] [ Modules ]

NAME

 defs_datatypes

FUNCTION

 This module contains definitions important structured datatypes for the ABINIT package.
 If you want to add one new datatype, please, examine first whether
 another datatype might meet your need (e.g. adding some records to it).
 Then, if you are sure your new structured datatype is needed,
 write it here, and DOCUMENT it properly (not all datastructure here are
 well documented, it is a shame ...).
 Do not forget: you will likely be the major winner if you document properly.
 Proper documentation of a structured datatype means:

  (1) Mention it in the list just below
  (2) Describe it in the NOTES section
  (3) Put it in alphabetical order in the the main section of this module
  (4) Document each of its records, except if they are described elsewhere
      (this exception is typically the case of the dataset associated with
      input variables, for which there is a help file)

 List of datatypes:
 * ebands_t: different information about the band structure
 * pseudopotential_type: for norm-conserving pseudopotential, all the information
 * pspheader_type: for norm-conserving pseudopotentials, the header of the file

COPYRIGHT

 Copyright (C) 2001-2024 ABINIT group (XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

35 #if defined HAVE_CONFIG_H
36 #include "config.h"
37 #endif
38 
39 #include "abi_common.h"
40 
41 module defs_datatypes
42 
43  use defs_basis
44 
45  implicit none

defs_datatypes/ebands_t [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 ebands_t

FUNCTION

 It contains different information about the band structure: eigenenergies, residuals, derivative of occupation
 numbers vs energy in case of metallic occupations and Brillouin zone according to the context: k points,
 occupation numbers, storage mode of wavefunctions, weights ...
 For example, the initial Brillouin zone, set up in the dataset, will be treated in the response function part of
 the code, to give a reduced Brillouin zone different from the original one, due to the breaking of the symmetries
 related to the existence of a wavevector, or the lack of time-reversal invariance.

SOURCE

 64  type ebands_t
 65 
 66   integer :: bantot                ! Total number of bands (sum(nband(:))
 67   integer :: ivalence              ! Highest valence band index (useful when occopt=9 only)
 68   integer :: mband                 ! Max number of bands i.e MAXVAL(nband) (to dimension arrays)
 69   integer :: nkpt                  ! Number of k points
 70   integer :: nspinor               ! 1 for collinear, 2 for noncollinear.
 71   integer :: nsppol                ! Number of spin-polarizations
 72   integer :: ntemp                 ! Number of temperatures
 73   integer :: occopt                ! Occupation option, see input variable.
 74 
 75   real(dp) :: entropy              ! Entropy associated with the smearing (adimensional)
 76   real(dp) :: fermie               ! Fermi energy ! CP: when occopt = 9, fermi energy of the quasi-FD distribution of excited
 77 ! electrons in the conduction bands above ivalence
 78   real(dp) :: fermih               ! Fermi energy of the excited holes in the valence bands <= ivalence (occopt = 9 only)
 79   real(dp) :: nelect               ! Number of electrons.
 80   real(dp) :: ne_qFD               ! Number of electrons excited in the bands > ivalence (occopt = 9 only)
 81   real(dp) :: nh_qFD               ! Number of holes     excited in the bands <=ivalence (occopt = 9 only)
 82   real(dp) :: tphysel              ! Physical temperature of electrons.
 83   real(dp) :: tsmear               ! Temperature of smearing.
 84   !real(dp) :: max_occ             ! Spin degeneracy factor: max_occ = two / (self%nspinor * self%nsppol)
 85 
 86   !real(dp) :: spinmagntarget
 87   ! TODO This should be set via dtset%spinmagntarget to simplify the API.
 88 
 89   integer,allocatable :: istwfk(:)
 90   ! istwfk(nkpt)
 91   ! Storage mode at each k point.
 92 
 93   integer,allocatable :: nband(:)
 94   ! nband(nkpt*nsppol)
 95   ! Number of bands at each k point and spin-polarisation.
 96 
 97   integer,allocatable :: npwarr(:)
 98   ! npwarr(nkpt)
 99   ! Number of plane waves at each k point.
100 
101   real(dp),allocatable :: kptns(:,:)
102   ! kptns(3,nkpt)
103   ! k-point vectors.
104 
105   real(dp),allocatable :: eig(:,:,:)
106   ! eig(mband,nkpt,nsppol)
107   ! Eigenvalues of each band.
108 
109   real(dp),allocatable :: linewidth(:,:,:,:)
110   ! linewidth(itemp,mband,nkpt,nsppol)
111   ! Linewidth of each band
112   ! MG: TODO: This array should be removed (I think Yannick introduced it, see also Ktmesh)
113 
114   real(dp),allocatable :: occ(:,:,:)
115   ! occ(mband, nkpt, nsppol)
116   ! occupation of each band.
117 
118   real(dp),allocatable :: doccde(:,:,:)
119   ! doccde(mband, nkpt, nsppol)
120   ! derivative of the occupation of each band wrt energy (needed for RF).
121 
122   real(dp),allocatable :: wtk(:)
123   ! wtk(nkpt)
124   ! weight of each k point, normalized to one.
125 
126   integer :: kptopt
127   ! Option used for k-point generation.
128 
129   integer :: nshiftk_orig, nshiftk
130   ! original number of shifts given in input and the actual value (changed in inkpts)
131 
132   real(dp) :: cellcharge
133   ! nelect = zion - cellcharge
134   ! Extra charge added to the unit cell when performing GS calculations
135   ! To treat a system missing one electron per unit cell, charge is set to +1.
136   ! When reading the band structure from an external file,
137   ! charge is mainly used as metadata describing the GS calculation that procuded the ebands_t object.
138   ! To simulate doping in a post-processing tool, use ebands_set_extrael that defines the value of %extra_el.
139   ! and changes %nelect, accordingly.
140 
141   real(dp) :: extrael = zero
142   ! Extra number of electrons.
143   ! This variable is mainly used to simulate doping in the rigid band approximation.
144   ! Set by ebands_set_extrael method.
145 
146   integer :: kptrlatt_orig(3,3), kptrlatt(3,3)
147   ! Original value of kptrlatt and value after the call to inkpts
148 
149   real(dp),allocatable :: shiftk_orig(:,:)
150   ! shiftk_orig(3, nshiftk_orig)
151   ! original shifts given in input (changed in inkpts).
152 
153   real(dp),allocatable :: shiftk(:,:)
154   ! shiftk(3, nshiftk)
155 
156  end type ebands_t

defs_datatypes/nctab_t [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 nctab_type

FUNCTION

  This object contains TABulated data for NC pseudopotentials
  It follows the conventions used in pawtab so that we can reuse
  the PAW routines for the treatment of model core change as well
  as the code in the atm2fft routines used to build approximated densities
  from atomic quantities. Methods are defined in m_psps.

SOURCE

222  type,public :: nctab_t
223 
224    integer :: mqgrid_vl = 0
225    ! Number of points in the reciprocal space grid on which
226    ! the radial functions are specified (same grid as the one used for the local part).
227 
228    ! TODO
229    !integer :: mqgrid_ff = 0
230 
231    logical :: has_tvale = .False.
232     ! True if the norm-conserving pseudopotential provides the atomic pseudized valence density.
233     ! If alchemy, has_tvale is True only if all the mixed pseudos
234     ! have the valence charge in the pseudopotential file.
235 
236    logical :: has_tcore = .False.
237     ! True if the norm-conserving pseudopotential has the model core-charge for NLCC.
238     ! If alchemy, has_tcore is True if at least one of the mixed pseudos has NLCC.
239     ! See also tcorespl
240 
241    real(dp) :: dncdq0 = zero
242     ! Gives 1/q d(tNcore(q))/dq for q=0
243     ! (tNcore(q) = FT of pseudo core density)
244 
245    real(dp) :: d2ncdq0 = zero
246     ! Gives contribution of d2(tNcore(q))/d2q for q=0
247     ! \int{(16/15)*pi^5*n(r)*r^6* dr}
248     ! (tNcore(q) = FT of pseudo core density)
249 
250    real(dp) :: dnvdq0 = zero
251     ! Gives 1/q d(tNvale(q))/dq for q=0
252     ! (tNvale(q) = FT of pseudo valence density)
253 
254    real(dp), allocatable :: tvalespl(:,:)
255     ! tvalespl(mqgrid_vl,2)
256     ! Gives the pseudo valence density in reciprocal space on a regular grid
257     ! Initialized only if has_nctvale(itypat)
258 
259    real(dp), allocatable :: tcorespl(:,:)
260     ! tcorespl(mqgrid_vl,2)
261     ! Gives the pseudo core density in reciprocal space on a regular grid.
262     ! tcorespl is **always** allocated and initialized with zeros if not has_tcore
263     ! A similar approach is used in PAW.
264 
265    integer :: num_tphi = 0
266    ! Number of pseudo atomic orbitals. 0 if pseudo does not provide them
267 
268    logical :: has_jtot = .False.
269    ! True if tpsi are given in terms of j (relativistic pseudo with SOC)
270 
271    real(dp), allocatable :: tphi_qspl(:,:,:)
272     ! (mqgrid_ff, 2, num_tphi)
273     ! Form factors for the pseudo wavefunctions.
274 
275    integer,allocatable :: tphi_n(:), tphi_l(:)
276     ! (num_tphi) arrays giving n, l
277 
278    real(dp),allocatable :: tphi_jtot(:)
279     ! (num_tphi) array with jtot.
280 
281    real(dp),allocatable :: tphi_occ(:)
282     ! (num_tphi) array with atomic occupancies taken from pseudo.
283 
284  end type nctab_t

defs_datatypes/pseudopotential_gth_type [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 pseudopotential_gth_type

FUNCTION

 This structure is a sub-structure of pseudopotential_type used to
 store parameters from the GTH pseudo-potentials. All arrays have
 indices running on 1:npsp for each read pseudo-file. The 'set' array
 is a check array, since several different pseudo can be used in a simulation
 it set a flag for each npsp if params have been set or not. This is
 redundant with psps%pspcod in the way that when psps%pspcod(i) is 2,
 then gth_params%set(i) is .true.. GTH pseudo previous to wavelets introduction
 doesn't have geometric information. These have been added on the last line.
 It is three radius information, the %hasGeometry flag is there to know
 which kind of pseudo has been read.

SOURCE

179  type pseudopotential_gth_type
180 
181 ! WARNING: if you modify this datatype, please check whether there might be creation/destruction/copy routines,
182 ! declared in another part of ABINIT, that might need to take into account your modification.
183 
184   real(dp), allocatable :: psppar(:, :, :)
185    ! These are {rloc, C(1...4)} coefficients for psppar(0, :, :) indices,
186    ! Followed by the h coefficients for psppar(1:2, :, :) indices.
187    !  size (0:2, 0:4, npsp)
188 
189   real(dp), allocatable :: radii_cf(:, :)
190    ! Cut-off radii for core part and long-range part.
191    ! radii_cf(:, 1) is for the long-range cut-off and
192    ! radii_cf(:, 2) is for the core cut-off. size (npsp, 2)
193 
194   real(dp), allocatable :: psp_k_par(:, :, :)
195    ! Spin orbit coefficients in HGH/GTH formats: k11p etc... see psp3ini.F90
196    !   dimension = num l channels, 3 coeffs, num psp = (1:lmax+1,1:3,npsp)
197 
198   logical, allocatable :: hasGeometry(:)
199    ! Flag for geometric information in the pseudo. size (npsp)
200 
201   logical, allocatable :: set(:)
202    ! Consistency array, used for checking size (npsp)
203 
204  end type pseudopotential_gth_type

defs_datatypes/pseudopotential_type [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 pseudopotential_type

FUNCTION

 This structured datatype contains all the information about one
 norm-conserving pseudopotential, including the description of the local
 and non-local parts, the different projectors, the non-linear core correction ...

SOURCE

300  type pseudopotential_type
301 
302 ! WARNING: if you modify this datatype, please check whether there might be creation/destruction/copy routines,
303 ! declared in another part of ABINIT, that might need to take into account your modification.
304 
305 ! Integer scalars
306   integer :: dimekb
307    ! Dimension of Ekb
308    ! ->Norm conserving : Max. number of Kleinman-Bylander energies for each atom type
309    !                     dimekb=lnmax (lnmax: see this file)
310    ! ->PAW : Max. number of Dij coefficients connecting projectors for each atom type
311    !                     dimekb=lmnmax*(lmnmax+1)/2 (lmnmax: see this file)
312 
313   integer :: lmnmax
314    !  If useylm=0, max number of (l,n)   comp. over all type of psps (lnproj)
315    !  If useylm=1, max number of (l,m,n) comp. over all type of psps (lmnproj)
316    !  If mpspso is 2, lmnmax takes into account the spin-orbit projectors,
317    !  so, it is equal to the max of lmnprojso or lnprojso, see pspheader_type
318 
319   integer :: lnmax
320    !  Max. number of (l,n) components over all type of psps
321    !  If mpspso is 2, lmnmax takes into account the spin-orbit projectors,
322    !  so, it is equal to the max of lnprojso, see pspheader_type
323 
324   integer :: mproj
325    ! Maximum number of non-local projectors over all angular momenta and type of psps
326    ! 0 only if all psps are local
327 
328   integer :: mpsang
329    ! Highest angular momentum of non-local projectors over all type of psps.
330    ! shifted by 1: for all local psps, mpsang=0; for largest s, mpsang=1,
331    ! for largest p, mpsang=2; for largest d, mpsang=3; for largest f, mpsang=4
332    ! This gives also the number of non-local "channels"
333 
334   integer :: mpspso
335    ! mpspso is set to 1 if none of the psps is used with a spin-orbit part (that
336    !  is, if the user input variable so_psp is not equal to 1 in at least one case
337    ! otherwise, it is set to 2
338 
339   integer :: mpssoang
340    ! Maximum number of channels, including those for treating the spin-orbit coupling
341    ! when mpspso=1, mpssoang=mpsang
342    ! when mpspso=2, mpssoang=2*mpsang-1
343 
344   integer :: mqgrid_ff
345    ! Number of points in the reciprocal space grid on which
346    ! the radial functions ffspl are specified
347 
348   integer :: mqgrid_vl
349    ! Number of points in the reciprocal space grid on which
350    ! the radial functions vlspl are specified
351 
352   integer :: mtypalch
353    ! Maximum number of alchemical pseudo atoms. If non-zero,
354    ! the mechanism to generate mixing of pseudopotentials is activated
355 
356   integer :: npsp
357    ! Number of types of pseudopotentials
358 
359   integer :: npspalch
360    ! Number of types of pseudopotentials used for alchemical purposes
361 
362   integer :: ntypat
363    ! Number of types of atoms (might be alchemy wrt pseudopotentials)
364 
365   integer :: ntypalch
366    ! Number of types of alchemical pseudoatoms
367 
368   integer :: ntyppure
369    ! Number of types of pure pseudoatoms
370 
371   integer :: n1xccc
372    ! Number of radial points for the description of the pseudo-core charge
373    ! (in the framework of the non-linear XC core correction)
374 
375   integer :: optnlxccc
376    ! Option for the choice of non-linear XC core correction treatment (see the input variable)
377    ! Used only for FHI pseudos.
378 
379   integer :: positron
380    ! Option for the choice of type of GS calculation (electron or positron)
381 
382   integer :: usepaw
383    ! if usepaw=0 , use norm-conserving psps part of the code
384    ! is usepaw=1 , use paw part of the code
385 
386   integer :: usewvl
387    ! if usewvl=0 ,  plane waves
388    ! is usewvl=1 ,  wavelets
389 
390   integer :: useylm
391    ! governs the way the nonlocal operator is to be applied:
392    !   1=using Ylm, 0=using Legendre polynomials
393 
394 ! Logical scalars
395   logical :: vlspl_recipSpace
396    ! governs if vlspl is compute in reciprocal space or in real
397    ! space (when available).
398 
399 ! Integer arrays
400   integer, allocatable :: algalch(:)
401    ! algalch(ntypalch)
402    ! For each type of pseudo atom, the algorithm to mix the pseudopotentials
403 
404   integer, allocatable :: indlmn(:,:,:)
405    ! indlmn(6,lmnmax,ntypat)
406    ! For each type of psp,
407    ! array giving l,m,n,lm,ln,spin for i=ln  (if useylm=0)
408    !                                or i=lmn (if useylm=1)
409    ! NB: spin is used for NC pseudos with SOC term: 1 if scalar term (spin diagonal), 2 if SOC term.
410 
411   integer, allocatable :: pspdat(:)
412    ! pspdat(ntypat)
413    ! For each type of psp, the date of psp generation, as given by the psp file
414 
415   integer, allocatable :: pspcod(:)
416    ! pspcod(npsp)
417    ! For each type of psp, the format -or code- of psp generation,
418    !  as given by the psp file
419 
420   integer, allocatable :: pspso(:)
421    ! pspso(npsps)
422    ! For each type of psp,
423    ! 1 if no spin-orbit component is taken into account
424    ! 2 if a spin-orbit component is used
425 
426   integer, allocatable :: pspxc(:)
427    ! pspxc(npsps)
428    ! For each type of psp, the XC functional that was used to generate it,
429    ! as given by the psp file
430 
431 ! Real (real(dp)) arrays
432 
433   real(dp), allocatable :: ekb(:,:)
434    ! ekb(dimekb,ntypat*(1-usepaw))
435    ! NORM-CONSERVING PSPS ONLY:
436    !    (Real) Kleinman-Bylander energies (hartree)
437    !           for number of basis functions (l,n) (lnmax)
438    !           and number of atom types (ntypat)
439    ! NOTE (MT):
440    !   ekb (norm-conserving) is now diagonal (one dimension lnmax);
441    !   it would be easy to give it a second (symmetric) dimension by putting
442    !   dimekb=lnmax*(lnmax+1)/2 in the place of dimekb=lmnmax.
443 
444   real(dp), allocatable :: ffspl(:,:,:,:)
445    ! ffspl(mqgrid_ff,2,lnmax,ntypat)
446    ! Gives, on the radial grid, the different non-local projectors,
447    ! in both the norm-conserving case, and the PAW case
448 
449   real(dp), allocatable :: mixalch(:,:)
450    ! mixalch(npspalch,ntypalch)
451    ! Mixing coefficients to generate alchemical pseudo atoms
452 
453   real(dp), allocatable :: qgrid_ff(:)
454    ! qgrid_ff(mqgrid_ff)
455    ! The coordinates of all the points of the radial grid for the nl form factors
456 
457   real(dp), allocatable :: qgrid_vl(:)
458    ! qgrid_vl(mqgrid_vl)
459    ! The coordinates of all the points of the radial grid for the local part of psp
460 
461   real(dp), allocatable :: vlspl(:,:,:)
462    ! vlspl(mqgrid_vl,2,ntypat)
463    ! Gives, on the radial grid, the local part of each type of psp.
464 
465   real(dp), allocatable :: dvlspl(:,:,:)
466    ! dvlspl(mqgrid_vl,2,ntypat)
467    ! Gives, on the radial grid, the first derivative of the local
468    ! part of each type of psp (computed when the flag 'vlspl_recipSpace' is true).
469 
470   real(dp), allocatable :: xcccrc(:)
471    ! xcccrc(ntypat)
472    ! Gives the maximum radius of the pseudo-core charge, for each type of psp.
473 
474   real(dp), allocatable :: xccc1d(:,:,:)
475    ! xccc1d(n1xccc*(1-usepaw),6,ntypat)
476    ! Norm-conserving psps only
477    ! The component xccc1d(n1xccc,1,ntypat) is the pseudo-core charge
478    ! for each type of atom, on the radial grid. The components
479    ! xccc1d(n1xccc,ideriv,ntypat) give the ideriv-th derivative of the
480    ! pseudo-core charge with respect to the radial distance.
481 
482   real(dp), allocatable :: zionpsp(:)
483    ! zionpsp(npsp)
484    ! For each pseudopotential, the ionic pseudo-charge
485    ! (giving raise to a long-range coulomb potential)
486 
487   real(dp), allocatable :: ziontypat(:)
488    ! ziontypat(ntypat)
489    ! For each type of atom (might be alchemy wrt psps), the ionic pseudo-charge
490    ! (giving raise to a long-range coulomb potential)
491 
492   real(dp), allocatable :: znuclpsp(:)
493    ! znuclpsp(npsp)
494    ! The atomic number of each pseudopotential
495 
496   real(dp), allocatable :: znucltypat(:)
497    ! znucltypat(ntypat)
498    ! The atomic number of each type of atom (might be alchemy wrt psps)
499 
500 ! Character arrays
501   character(len=fnlen), allocatable :: filpsp(:)
502    ! filpsp(npsps)
503    ! The filename of the pseudopotential
504 
505   character(len=fnlen), allocatable :: title(:)
506    ! title(npsp)
507    ! The content of first line read from the psp file
508 
509   character(len=md5_slen), allocatable :: md5_pseudos(:)
510   ! md5pseudos(npsp)
511   ! md5 checksums associated to pseudos (read from file)
512 
513   type(pseudopotential_gth_type) :: gth_params
514    ! Types for pseudo-potentials that are based on parameters. Currently, only
515    ! GTH are supported (see pseudopotential_gth_type). To add one, one should
516    ! create an initialisation method and a destruction method in 02psp (see
517    ! psp2params.F90). These methods are called in driver().
518 
519    type(nctab_t),allocatable :: nctab(:)
520    ! nctab(ntypat)
521    ! Tables storing additional data for NC pseudopotentials that are not always avaiable if every psp format.
522    ! We try to mimim pawtab as much as possible so that we can reuse PAW routines in the NC context.
523 
524    integer :: nc_xccc_gspace = 0
525    ! NC pseudos only. Set to 1 if the non-linear core correction should
526    ! be treated in G-space similarly to the approach used for PAW.
527 
528  end type pseudopotential_type

defs_datatypes/pspheader_paw_type [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 pspheader_paw_type

FUNCTION

 The pspheader_paw_type structured datatype gather additional information
 about a PAW pseudopotential file, from its header.

SOURCE

543  type pspheader_paw_type
544 
545 ! WARNING: if you modify this datatype, please check whether there might be creation/destruction/copy routines,
546 ! declared in another part of ABINIT, that might need to take into account your modification.
547 
548 ! WARNING: Also pay attention to subroutine pspheads_comm, which broadcasts this datatype.
549 
550   integer :: basis_size    ! Number of elements of the wf basis ((l,n) quantum numbers)
551   integer :: l_size        ! Maximum value of l+1 leading to a non zero Gaunt coefficient
552   integer :: lmn_size      ! Number of elements of the paw basis
553   integer :: mesh_size     ! Dimension of (main) radial mesh
554   integer :: pawver        ! Version number of paw psp format
555   integer :: shape_type    ! Type of shape function
556   real(dp) :: rpaw         ! Radius for paw spheres
557   real(dp) :: rshp         ! Cut-off radius of shape function
558 
559  end type pspheader_paw_type

defs_datatypes/pspheader_type [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 pspheader_type

FUNCTION

 The pspheader_type structured datatype gather different information
 about a pseudopotential file, from its header.

SOURCE

574  type pspheader_type
575 
576 ! WARNING: if you modify this datatype, please check whether there might be creation/destruction/copy routines,
577 ! declared in another part of ABINIT, that might need to take into account your modification.
578 ! WARNING: Also pay attention to subroutine pspheads_comm, which broadcasts this datatype.
579 
580   integer :: nproj(0:3) ! number of scalar projectors for each angular momentum
581 
582   integer :: nprojso(3) ! number of spin-orbit projectors for each angular momentum
583 
584   integer :: lmax       ! maximum l quantum number (-1 if only local)
585                         ! Example : s only       -> lmax=0
586                         !           s and p      -> lmax=1
587                         !           s,p,d        -> lmax=2
588 
589   integer :: pspcod
590     ! code number of the pseudopotential
591 
592   integer :: pspdat
593     ! date of generation of the pseudopotential
594 
595   integer :: pspxc
596     ! exchange-correlation functional
597 
598   integer :: pspso
599     ! spin-orbit characteristics
600     ! 0 if pseudo does not support SOC, 1 or 2 if SOC terms are available in the pp file.
601     ! Note that one could have a pseudo with SOC terms but ignore the SOC contribution
602     ! For example, it's possible to use nspinor=2 and set so_psp to 0 in the input file
603     ! or perform a run with nspinor=1 with pseudos containing SOC terms.
604 
605   integer :: usewvl
606    ! if usewvl=0 ,  plane waves
607    ! is usewvl=1 ,  wavelets
608 
609   integer :: xccc
610     ! =0 if no XC core correction, non-zero if XC core correction
611 
612   real(dp) :: zionpsp
613     ! charge of the ion made of core electrons only
614 
615   real(dp) :: znuclpsp
616     ! atomic number of the nuclei
617 
618   real(dp) :: GTHradii(0:4)
619     ! Radii values for GTH (and HGH) family potentials
620 
621   character(len=fnlen) :: filpsp
622     ! name of the psp file
623 
624   character(len=fnlen) :: title
625     ! content of first line read from the psp file
626 
627   character(len=md5_slen) :: md5_checksum = md5_none
628     ! md5 checksum read from file.
629 
630   type(pspheader_paw_type) :: pawheader
631     ! only for PAW psps ; see m_pawpsp.
632 
633  end type pspheader_type