TABLE OF CONTENTS


ABINIT/m_psps [ Modules ]

[ Top ] [ Modules ]

NAME

  m_psps

FUNCTION

  This module provides method to allocate/free/initialize the
  pseudopotential_type object.

COPYRIGHT

  Copyright (C) 2014-2024 ABINIT group (XG,DC,MG)
  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

17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 module m_psps
24 
25  use defs_basis
26  use m_abicore
27  use m_errors
28  use m_xmpi
29  use m_nctk
30  use m_copy
31  use m_dtset
32 #ifdef HAVE_NETCDF
33  use netcdf
34 #endif
35 
36  use m_fstrings,      only : itoa, sjoin, yesno, atoi
37  use m_io_tools,      only : open_file
38  use m_symtk,         only : matr3inv
39  use defs_datatypes,  only : pspheader_type, pseudopotential_type, pseudopotential_gth_type, nctab_t
40  use m_paw_numeric,   only : paw_spline
41  use m_pawrad,        only : pawrad_type, pawrad_init, pawrad_free, simp_gen
42  use m_pawpsp,        only : pawpsp_cg
43  use m_parser,        only : chkint_eq
44  use m_memeval,       only : getdim_nloc, setmqgrid
45 
46  implicit none
47 
48  private
49 
50  ! Helper functions
51  public :: test_xml_xmlpaw_upf     ! Test if a pseudo potential file is in XML, XML-PAW or in UPF format.
52 
53  public :: psps_init_global        ! Allocate and init all part of psps structure that are independent of a given dataset.
54  public :: psps_init_from_dtset    ! Allocate and init all part of psps structure that are dependent of a given dataset.
55  public :: psps_free               ! Deallocate all memory of psps structure.
56  public :: psps_copy               ! Copy the psps structure.
57  public :: psps_print              ! Print info on the pseudopotentials.
58  public :: psps_ncwrite_path       ! Create a netcdf file and write psps data.
59  public :: psps_ncwrite            ! Write psps data in an open netcdf file.
60  public :: psps_ncread             ! Read psps data from an open netcdf file.
61 
62  public :: nctab_init              ! Create the object.
63  public :: nctab_free              ! Free memory.
64  public :: nctab_copy              ! Copy the object.
65  public :: nctab_eval_tvalespl     ! Evaluate spline-fit of the atomic pseudo valence charge in reciprocal space.
66  public :: nctab_eval_tcorespl     ! Evalute spline-fit of the model core charge in reciprocal space.
67  public :: nctab_mixalch           ! Mix the pseudopotential tables. Used for alchemical mixing.

m_psps/nctab_copy [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  nctab_copy

FUNCTION

  Copy the object.

SOURCE

1491 subroutine nctab_copy(nctabin, nctabout)
1492 
1493 !Arguments ------------------------------------
1494  class(nctab_t),intent(in) :: nctabin
1495  class(nctab_t),intent(inout) :: nctabout
1496 
1497 ! *************************************************************************
1498 
1499  nctabout%mqgrid_vl  = nctabin%mqgrid_vl
1500  nctabout%has_tvale  = nctabin%has_tvale
1501  nctabout%has_tcore  = nctabin%has_tcore
1502  nctabout%dncdq0     = nctabin%dncdq0
1503  nctabout%d2ncdq0    = nctabin%d2ncdq0
1504  nctabout%dnvdq0     = nctabin%dnvdq0
1505 
1506  ! TODO Why not check values of has_tvale and has_tcore?
1507  if (allocated(nctabin%tvalespl)) call alloc_copy(nctabin%tvalespl, nctabout%tvalespl)
1508  if (allocated(nctabin%tcorespl)) call alloc_copy(nctabin%tcorespl, nctabout%tcorespl)
1509 
1510 end subroutine nctab_copy

m_psps/nctab_eval_tcorespl [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  nctab_eval_tcorespl

FUNCTION

  Evalute spline-fit of the model core charge in reciprocal space.

INPUTS

  xcccrc=maximum radius of the pseudo-core charge
  n1xccc=Number of radial points for the description of the pseudo-core charge
     (in the framework of the non-linear XC core correction)
  mqgrid_vl=Number of points in the reciprocal space grid
  qgrid_vl(mqgrid_vl)=The coordinates of all the points of the radial q-grid
  xccc1d(n1xccc,6)= The component xccc1d(n1xccc,1) is the pseudo-core charge
   on the radial grid. The components xccc1d(n1xccc,ideriv) give the ideriv-th derivative of the
   pseudo-core charge with respect to the radial distance.

SIDE EFFECTS

  nctabl%tcorespl(mqgrid_vl,2)
  nctab%d2ncdq0
  nctab%dncdq0

SOURCE

1597 subroutine nctab_eval_tcorespl(nctab, n1xccc, xcccrc, xccc1d, mqgrid_vl, qgrid_vl)
1598 
1599 !Arguments ------------------------------------
1600 !scalars
1601  class(nctab_t),intent(inout) :: nctab
1602  integer,intent(in) :: n1xccc,mqgrid_vl
1603  real(dp),intent(in) :: xcccrc
1604 !arrays
1605  real(dp),intent(in) :: xccc1d(n1xccc,6),qgrid_vl(mqgrid_vl)
1606 
1607 !Local variables-------------------------------
1608  real(dp) :: amesh,yp1,ypn
1609  type(pawrad_type) :: core_mesh
1610 
1611 ! *************************************************************************
1612 
1613  ABI_CHECK(mqgrid_vl == nctab%mqgrid_vl, "wrong mqgrid_vl")
1614 
1615  if (.not. allocated(nctab%tcorespl)) then
1616    ABI_CALLOC(nctab%tcorespl, (mqgrid_vl, 2))
1617  else
1618    ABI_CHECK(size(nctab%tcorespl, dim=1) == mqgrid_vl, "wrong mqgrid_vl")
1619  end if
1620 
1621  ! Skip loop if this atom has no core charge
1622  if (abs(xcccrc) < tol16) then
1623    nctab%has_tcore = .False.
1624    return
1625  end if
1626 
1627  nctab%has_tcore = .True.
1628  ! XCCC is given on a linear mesh.
1629  amesh = xcccrc / dble(n1xccc-1)
1630  call pawrad_init(core_mesh, mesh_size=n1xccc, mesh_type=1, rstep=amesh)
1631 
1632  ! Compute 4\pi\int[(\frac{\sin(2\pi q r)}{2\pi q r})(r^2 n(r))dr].
1633  ! write(std_out,*)"xccc1d: amesh, min, max, minloc ",amesh,maxval(xccc1d(:,1)),minval(xccc1d(:,1)),minloc(xccc1d(:,1))
1634  call pawpsp_cg(nctab%dncdq0, nctab%d2ncdq0, mqgrid_vl, qgrid_vl, nctab%tcorespl(:,1), &
1635                 core_mesh, xccc1d(:,1), yp1, ypn)
1636 
1637  ! Compute second derivative of tcorespl(q)
1638  call paw_spline(qgrid_vl, nctab%tcorespl(:,1), mqgrid_vl, yp1, ypn, nctab%tcorespl(:,2))
1639  call pawrad_free(core_mesh)
1640 
1641 end subroutine nctab_eval_tcorespl

m_psps/nctab_eval_tvalespl [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  nctab_eval_tvalespl

FUNCTION

  Evalute spline-fit of the atomic pseudo valence charge in reciprocal space.

INPUTS

  zion=nominal valence of atom as specified in psp file. Used to rescale the f(q=0) component
  mesh<pawrad_type>Radial mesh (r-space) used for the valence denity.
  valr(mesh%mesh_size)=Valence density in real space.
  mqgrid_vl=Number of points in the reciprocal space grid
  qgrid_vl(mqgrid_vl)=The coordinates of all the points of the radial q-grid

SIDE EFFECTS

  nctabl%tvalspl(mqgrid_vl,2)
  nctab%d2ncdq0

SOURCE

1533 subroutine nctab_eval_tvalespl(nctab, zion, mesh, valr, mqgrid_vl, qgrid_vl)
1534 
1535 !Arguments ------------------------------------
1536  class(nctab_t),intent(inout) :: nctab
1537  integer,intent(in) :: mqgrid_vl
1538  real(dp),intent(in) :: zion
1539  type(pawrad_type),intent(in) :: mesh
1540 !arrays
1541  real(dp),intent(in) :: valr(mesh%mesh_size),qgrid_vl(mqgrid_vl)
1542 
1543 !Local variables-------------------------------
1544  real(dp) :: fact,yp1,ypn,d2nvdq0
1545 
1546 ! *************************************************************************
1547 
1548  nctab%has_tvale = .True.
1549  if (.not. allocated(nctab%tvalespl)) then
1550    ABI_MALLOC(nctab%tvalespl, (mqgrid_vl, 2))
1551  else
1552    ABI_CHECK(size(nctab%tvalespl, dim=1) == mqgrid_vl, "wrong mqgrid_vl")
1553  end if
1554 
1555  call pawpsp_cg(nctab%dnvdq0, d2nvdq0, mqgrid_vl, qgrid_vl, nctab%tvalespl(:,1), mesh, valr, yp1, ypn)
1556  call simp_gen(yp1, mesh%rad**2 * valr, mesh)
1557  write(std_out,*)" valence charge (before rescaling) integrates to: ",four_pi*yp1
1558 
1559  ! Rescale the integral to have the correct number of valence electrons.
1560  ! In some cases, indeed, the radial mesh is not large enough and some valence charge is missing
1561  ! pawpsp_cg extrapolates the integrand beyond rmax but this is not enough.
1562  ! Remember that tvalespl is used to build an initial guess for rhor hence it's very important
1563  ! to have the correct electrostatic.
1564  fact = zion / nctab%tvalespl(1,1)
1565  nctab%tvalespl(:,1) = nctab%tvalespl(:,1) * fact
1566 
1567  ! Compute second derivative of tvalespl(q)
1568  call paw_spline(qgrid_vl,nctab%tvalespl(:,1),mqgrid_vl,yp1,ypn,nctab%tvalespl(:,2))
1569 
1570 end subroutine nctab_eval_tvalespl

m_psps/nctab_free [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  nctab_free

FUNCTION

 Free memory allocated in nctab_t

SOURCE

1464 subroutine nctab_free(nctab)
1465 
1466 !Arguments ------------------------------------
1467  class(nctab_t),intent(inout) :: nctab
1468 
1469 ! *************************************************************************
1470 
1471  ABI_SFREE(nctab%tvalespl)
1472  ABI_SFREE(nctab%tcorespl)
1473  ABI_SFREE(nctab%tphi_qspl)
1474  ABI_SFREE(nctab%tphi_n)
1475  ABI_SFREE(nctab%tphi_l)
1476  ABI_SFREE(nctab%tphi_jtot)
1477  ABI_SFREE(nctab%tphi_occ)
1478 
1479 end subroutine nctab_free

m_psps/nctab_init [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  nctab_init

FUNCTION

  Create nctab_t.

INPUTS

  mqgrid_vl=Number of q-points
  has_tcore=True if the pseudo has NLCC.
  has_tvale=True if the atomic valence density is available.

SOURCE

1427 subroutine nctab_init(nctab, mqgrid_vl, has_tcore, has_tvale)
1428 
1429 !Arguments ------------------------------------
1430  class(nctab_t),intent(inout) :: nctab
1431  integer,intent(in) :: mqgrid_vl
1432  logical,intent(in) :: has_tcore, has_tvale
1433 
1434 ! *************************************************************************
1435 
1436  nctab%mqgrid_vl = mqgrid_vl
1437 
1438  ! The array for the model core charge is always allocated and initialized with zeros.
1439  ! This approach is similar to the one used in the PAW code.
1440  ! has_tcore tells us whether the model core charge is present or not.
1441  nctab%has_tcore = has_tcore
1442  nctab%dncdq0 = zero; nctab%d2ncdq0 = zero
1443  ABI_CALLOC(nctab%tcorespl, (mqgrid_vl, 2))
1444 
1445  ! tvalespl is allocated only if available.
1446  nctab%has_tvale = has_tvale
1447  nctab%dnvdq0 = zero
1448  if (has_tvale) then
1449    ABI_CALLOC(nctab%tvalespl, (mqgrid_vl, 2))
1450  end if
1451 
1452 end subroutine nctab_init

m_psps/nctab_mixalch [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  nctab_mixalch

FUNCTION

  Mix the pseudopotential tables. Used for alchemical mixing.

INPUTS

 nctabs(npspalch)=NC tables to be mixed
 npspalch=Number of alchemical pseudos.
 ntypalch=Number of types of alchemical pseudoatoms
 algalch(ntypalch)=For each type of pseudo atom, the algorithm to mix the pseudopotentials
 mixalch(npspalch,ntypalch)=Mixing coefficients to generate alchemical pseudo atoms

OUTPUT

 mixtabs(ntypalch)=NC tables describing the alchemical pseudos

SOURCE

1663 subroutine nctab_mixalch(nctabs, npspalch, ntypalch, algalch, mixalch, mixtabs)
1664 
1665 !Arguments ------------------------------------
1666 !scalars
1667  integer,intent(in) :: npspalch,ntypalch
1668 !arrays
1669  integer,intent(in) :: algalch(ntypalch)
1670  real(dp),intent(in) :: mixalch(npspalch, ntypalch)
1671  type(nctab_t),intent(in) :: nctabs(npspalch)
1672  type(nctab_t),target,intent(inout) :: mixtabs(ntypalch)
1673 
1674 !Local variables-------------------------------
1675 !scalars
1676  integer :: ipspalch,itypalch
1677  logical :: has_tcore, has_tvale
1678  real(dp) :: mc
1679  type(nctab_t),pointer :: mix
1680 
1681 ! *************************************************************************
1682 
1683  ABI_CHECK(all(nctabs(:)%mqgrid_vl == nctabs(1)%mqgrid_vl), "Wrong mqgrid_vl")
1684  ABI_CHECK(all(algalch == 1), "algalch /= 1 not implemented")
1685 
1686  do itypalch=1,ntypalch
1687 
1688    ! has_tcore is true is at least one pseudo has nlcc.
1689    ! has_tvale is true if *all* mixed pseudos have the PS valence charge.
1690    has_tcore = .False.; has_tvale = .True.
1691    do ipspalch=1,npspalch
1692      if (abs(mixalch(ipspalch,itypalch)) < tol6) cycle
1693      if (nctabs(ipspalch)%has_tcore) has_tcore = .True.
1694      if (.not. nctabs(ipspalch)%has_tvale) has_tvale = .False.
1695    end do
1696    !write(std_out,*)has_tvale, has_tcore
1697 
1698    call nctab_free(mixtabs(itypalch))
1699    call nctab_init(mixtabs(itypalch), nctabs(1)%mqgrid_vl, has_tcore, has_tvale)
1700    mix => mixtabs(itypalch)
1701 
1702    do ipspalch=1,npspalch
1703      mc = mixalch(ipspalch,itypalch)
1704      if (abs(mc) < tol6) cycle
1705      ! Linear combination of the quantities
1706      ! Mix core for NLCC
1707      if (has_tcore) then
1708        mix%tcorespl = mix%tcorespl + mc * nctabs(ipspalch)%tcorespl
1709        mix%dncdq0 = mix%dncdq0 + mc * nctabs(ipspalch)%dncdq0
1710        mix%d2ncdq0 = mix%d2ncdq0 + mc * nctabs(ipspalch)%d2ncdq0
1711      end if
1712      ! Mix pseudo valence charge.
1713      if (has_tvale) then
1714        mix%tvalespl = mix%tvalespl + mc * nctabs(ipspalch)%tvalespl
1715        mix%dnvdq0 = mix%dnvdq0 + mc * nctabs(ipspalch)%dnvdq0
1716      end if
1717    end do
1718 
1719  end do
1720 
1721 end subroutine nctab_mixalch

m_psps/psp2params_copy [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psp2params_copy

FUNCTION

INPUTS

OUTPUT

SOURCE

1347 subroutine psp2params_copy(gth_paramsin, gth_paramsout)
1348 
1349 !Arguments ------------------------------------
1350  class(pseudopotential_gth_type),intent(in) :: gth_paramsin
1351  class(pseudopotential_gth_type),intent(inout) :: gth_paramsout
1352 
1353 ! *********************************************************************
1354 
1355  if (allocated(gth_paramsin%psppar)) then
1356    call alloc_copy( gth_paramsin%psppar, gth_paramsout%psppar)
1357  end if
1358  if (allocated(gth_paramsin%radii_cf)) then
1359    call alloc_copy( gth_paramsin%radii_cf, gth_paramsout%radii_cf)
1360  end if
1361  if (allocated(gth_paramsin%psp_k_par)) then
1362    call alloc_copy( gth_paramsin%psp_k_par, gth_paramsout%psp_k_par)
1363  end if
1364  if (allocated(gth_paramsin%hasGeometry)) then
1365    call alloc_copy( gth_paramsin%hasGeometry, gth_paramsout%hasGeometry)
1366  end if
1367  if (allocated(gth_paramsin%set)) then
1368    call alloc_copy( gth_paramsin%set, gth_paramsout%set)
1369  end if
1370 
1371 end subroutine psp2params_copy

m_psps/psp2params_free [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psp2params_free

FUNCTION

 Deallocate a previously allocated data structure for storage of GTH parameters.

INPUTS

SIDE EFFECTS

  gth_params <type (pseudopotential_gth_type)>=the values to deallocate.

SOURCE

1390 subroutine psp2params_free(gth_params)
1391 
1392 !Arguments ------------------------------------
1393 !scalars
1394  type(pseudopotential_gth_type),intent(inout) :: gth_params
1395 
1396 ! *********************************************************************
1397 
1398  ABI_SFREE(gth_params%set)
1399  ABI_SFREE(gth_params%hasGeometry)
1400 
1401  ! Coefficients for local part and projectors
1402  ABI_SFREE(gth_params%psppar)
1403 
1404  ! Coefficients for spin orbit part
1405  ABI_SFREE(gth_params%psp_k_par)
1406 
1407  ! Different radii
1408  ABI_SFREE(gth_params%radii_cf)
1409 
1410 end subroutine psp2params_free

m_psps/psp2params_init [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psp2params_init

FUNCTION

 Allocate and initialise the data structure holding parameters for the GTH
 pseudo-potentials.

  MJV note: this should be renamed: psp2 suggests it relates to pspcod 2,
     whereas it is actually 3
    the parameters would also be better off separated into C and h arrays

INPUTS

  npsp=number of true pseudo used (not alchemy).

OUTPUT

  gth_params <type (pseudopotential_gth_type)>=the values to allocate and initialise.

SOURCE

1302 subroutine psp2params_init(gth_params, npsp)
1303 
1304 !Arguments ------------------------------------
1305  class(pseudopotential_gth_type),intent(out) :: gth_params
1306  integer,intent(in) :: npsp
1307 
1308 ! *********************************************************************
1309 
1310 !Check array, no params are currently set.
1311  ABI_MALLOC(gth_params%set,(npsp))
1312  gth_params%set(:) = .false.
1313 
1314 !Check array, have geometric information been filled?
1315  ABI_MALLOC(gth_params%hasGeometry,(npsp))
1316  gth_params%hasGeometry(:) = .false.
1317 
1318 !Coefficients for local part and projectors
1319  ABI_MALLOC(gth_params%psppar,(0:4, 0:6, npsp))
1320  gth_params%psppar = zero
1321 
1322 !Coefficients for spin orbit part
1323  ABI_MALLOC(gth_params%psp_k_par,(1:4, 1:3, npsp))
1324  gth_params%psp_k_par = zero
1325 
1326 !Different radii
1327  ABI_MALLOC(gth_params%radii_cf,(npsp, 3))
1328  gth_params%radii_cf = zero
1329 
1330 end subroutine psp2params_init

m_psps/psps_copy [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_copy

FUNCTION

 Copy the psps structure.

SOURCE

549 subroutine psps_copy(pspsin, pspsout)
550 
551 !Arguments ------------------------------------
552  class(pseudopotential_type),intent(in) :: pspsin
553  class(pseudopotential_type),intent(inout) :: pspsout
554 
555 !Local variables-------------------------------
556  integer :: ii
557 
558 ! *************************************************************************
559 
560  ! integer
561  pspsout%dimekb         = pspsin%dimekb
562  pspsout%lmnmax         = pspsin%lmnmax
563  pspsout%lnmax          = pspsin%lnmax
564  pspsout%mproj          = pspsin%mproj
565  pspsout%mpsang         = pspsin%mpsang
566  pspsout%mpspso         = pspsin%mpspso
567  pspsout%mpssoang       = pspsin%mpssoang
568  pspsout%mqgrid_ff      = pspsin%mqgrid_ff
569  pspsout%mqgrid_vl      = pspsin%mqgrid_vl
570  pspsout%mtypalch       = pspsin%mtypalch
571  pspsout%npsp           = pspsin%npsp
572  pspsout%npspalch       = pspsin%npspalch
573  pspsout%ntypat         = pspsin%ntypat
574  pspsout%ntypalch       = pspsin%ntypalch
575  pspsout%ntyppure       = pspsin%ntyppure
576  pspsout%n1xccc         = pspsin%n1xccc
577  pspsout%optnlxccc      = pspsin%optnlxccc
578  pspsout%positron       = pspsin%positron
579  pspsout%usepaw         = pspsin%usepaw
580  pspsout%usewvl         = pspsin%usewvl
581  pspsout%useylm         = pspsin%useylm
582  pspsout%nc_xccc_gspace = pspsin%nc_xccc_gspace
583 
584  ! logical
585  pspsout%vlspl_recipSpace = pspsin%vlspl_recipSpace
586 
587  ! integer allocatable
588  if (allocated(pspsin%algalch)) call alloc_copy(pspsin%algalch, pspsout%algalch)
589  if (allocated(pspsin%indlmn)) call alloc_copy(pspsin%indlmn, pspsout%indlmn)
590  if (allocated(pspsin%pspdat)) call alloc_copy(pspsin%pspdat, pspsout%pspdat)
591  if (allocated(pspsin%pspcod)) call alloc_copy(pspsin%pspcod, pspsout%pspcod)
592  if (allocated(pspsin%pspso)) call alloc_copy(pspsin%pspso, pspsout%pspso)
593  if (allocated(pspsin%pspxc)) call alloc_copy(pspsin%pspxc, pspsout%pspxc)
594 
595  ! real allocatable
596  if (pspsin%dimekb > 0 .and. pspsin%usepaw==0) then
597    if (allocated(pspsin%ekb)) then
598      call alloc_copy( pspsin%ekb, pspsout%ekb)
599    end if
600  else
601    ABI_MALLOC(pspsout%ekb,(pspsout%dimekb,pspsout%ntypat * (1 - pspsout%usepaw)))
602    pspsout%ekb = zero
603  end if
604  if (allocated(pspsin%ffspl)) call alloc_copy( pspsin%ffspl, pspsout%ffspl)
605  if (allocated(pspsin%mixalch)) call alloc_copy(pspsin%mixalch, pspsout%mixalch)
606  if (allocated(pspsin%qgrid_ff)) call alloc_copy(pspsin%qgrid_ff, pspsout%qgrid_ff)
607  if (allocated(pspsin%qgrid_vl)) call alloc_copy(pspsin%qgrid_vl, pspsout%qgrid_vl)
608  if (allocated(pspsin%vlspl)) call alloc_copy(pspsin%vlspl, pspsout%vlspl)
609  if (allocated(pspsin%dvlspl)) call alloc_copy(pspsin%dvlspl, pspsout%dvlspl)
610 
611  if (allocated(pspsin%ziontypat)) call alloc_copy(pspsin%ziontypat, pspsout%ziontypat)
612  if (allocated(pspsin%znucltypat)) call alloc_copy(pspsin%znucltypat, pspsout%znucltypat)
613 
614  ! GA: Could make a check on mtypalch here
615  if (allocated(pspsin%znuclpsp)) call alloc_copy(pspsin%znuclpsp, pspsout%znuclpsp)
616  if (allocated(pspsin%zionpsp)) call alloc_copy(pspsin%zionpsp, pspsout%zionpsp)
617 
618  if (pspsin%n1xccc > 0) then
619    if (allocated(pspsin%xcccrc)) call alloc_copy(pspsin%xcccrc, pspsout%xcccrc)
620    if (allocated(pspsin%xccc1d)) call alloc_copy(pspsin%xccc1d, pspsout%xccc1d)
621  end if
622 
623  ! allocate and copy character strings
624  ABI_MALLOC(pspsout%filpsp,(pspsout%npsp))
625  ABI_MALLOC(pspsout%title,(pspsout%npsp))
626  ABI_MALLOC(pspsout%md5_pseudos,(pspsout%npsp))
627  do ii=1,pspsout%npsp
628    pspsout%filpsp(ii) = pspsin%filpsp(ii)
629    pspsout%title(ii) = pspsin%title(ii)
630    pspsout%md5_pseudos(ii) = pspsin%md5_pseudos(ii)
631  end do
632 
633  ! allocate and copy objects
634  if (allocated(pspsin%nctab)) then
635    ABI_MALLOC(pspsout%nctab,(pspsout%ntypat))
636    if (pspsin%usepaw==0) then
637      do ii=1,pspsout%ntypat
638        call nctab_copy(pspsin%nctab(ii), pspsout%nctab(ii))
639      end do
640    end if
641  end if
642 
643  call psp2params_copy(pspsin%gth_params, pspsout%gth_params)
644 
645 end subroutine psps_copy

m_psps/psps_free [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_free

FUNCTION

 Deallocate all memory of psps structure.

SOURCE

491 subroutine psps_free(psps)
492 
493 !Arguments ------------------------------------
494  type(pseudopotential_type),intent(inout) :: psps
495 
496 !Local variables-------------------------------
497  integer :: ii
498 
499 ! *************************************************************************
500 
501 !Allocation of some arrays independent of the dataset
502  ABI_SFREE(psps%filpsp)
503  ABI_SFREE(psps%pspcod)
504  ABI_SFREE(psps%pspdat)
505  ABI_SFREE(psps%pspso)
506  ABI_SFREE(psps%pspxc)
507  ABI_SFREE(psps%title)
508  ABI_SFREE(psps%algalch)
509  ABI_SFREE(psps%mixalch)
510  ABI_SFREE(psps%ekb)
511  ABI_SFREE(psps%indlmn)
512  ABI_SFREE(psps%ffspl)
513  ABI_SFREE(psps%qgrid_ff)
514  ABI_SFREE(psps%qgrid_vl)
515  ABI_SFREE(psps%vlspl)
516  ABI_SFREE(psps%dvlspl)
517  ABI_SFREE(psps%xccc1d)
518  ABI_SFREE(psps%xcccrc)
519  ABI_SFREE(psps%ziontypat)
520  ABI_SFREE(psps%zionpsp)
521  ABI_SFREE(psps%znucltypat)
522  ABI_SFREE(psps%znuclpsp)
523  ABI_SFREE(psps%md5_pseudos)
524 
525  ! Free types.
526  call psp2params_free(psps%gth_params)
527 
528  if (allocated(psps%nctab)) then
529    do ii=1,size(psps%nctab)
530      call nctab_free(psps%nctab(ii))
531    end do
532    ABI_FREE(psps%nctab)
533  end if
534 
535 end subroutine psps_free

m_psps/psps_init_from_dtset [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_init_from_dtset

FUNCTION

 Allocate and initialise all part of psps structure that are dependent of a given dataset.

INPUTS

 dtset=<type dataset_type>a given dataset
 pspheads(npsp)=<type pspheader_type>all the important information from the
   pseudopotential file header, as well as the psp file name

SIDE EFFECTS

 psps=<type pseudopotential_type>the pseudopotentials description

SOURCE

253 subroutine psps_init_from_dtset(psps, dtset, idtset, pspheads)
254 
255 !Arguments ------------------------------------
256 !scalars
257  class(pseudopotential_type),intent(inout) :: psps
258  integer,intent(in) :: idtset
259  type(dataset_type),intent(in) :: dtset
260 !arrays
261  type(pspheader_type),intent(in) :: pspheads(psps%npsp)
262 
263 !Local variables-------------------------------
264 !scalars
265  integer,save :: dimekb_old=-1,lmnmax_old=-1,lnmax_old=-1,mqgridff_old=0
266  integer,save :: mqgridvl_old=0,ntypat_old=-1,usepaw_old=-1
267  integer :: ipsp,lmnmax,lmnmaxso,lnmax,lnmaxso,newmqgrid,newmqgriddg,nptsgvec
268  integer :: changed,ii,itypat
269  real(dp) :: gprimd_orig(3,3)
270 
271 ! *************************************************************************
272 
273  psps%optnlxccc   = dtset%optnlxccc
274 !Determine the number of points needed in reciprocal space to represent the
275 !pseudopotentials (either set by hand from input variable or set automatically by abinit)
276  nptsgvec         = 200 !This has to be chosen one and for all or else ??
277  newmqgrid        = dtset%mqgrid
278  newmqgriddg      = dtset%mqgriddg
279 
280  !JB:Which image to use ? I guess 1 always works
281  call matr3inv(dtset%rprimd_orig(:,:,1),gprimd_orig)
282  if ( dtset%usewvl == 0) then
283    call setmqgrid(newmqgrid,newmqgriddg,dtset%ecut*dtset%dilatmx**2,&
284 &       dtset%pawecutdg*dtset%dilatmx**2,gprimd_orig,nptsgvec,psps%usepaw)
285  else
286    call setmqgrid(newmqgrid,newmqgriddg,one,one,gprimd_orig,nptsgvec,psps%usepaw)
287  end if
288  psps%mqgrid_ff   = newmqgrid
289  if (psps%usepaw == 1) then
290    psps%mqgrid_vl = newmqgriddg
291  else
292    psps%mqgrid_vl = newmqgrid
293  end if
294 
295 !Determine the maximum number of projectors, for the set of pseudo atom
296  call getdim_nloc(lmnmax,lmnmaxso,lnmax,lnmaxso,dtset%mixalch_orig,dtset%nimage,psps%npsp,dtset%npspalch,&
297 & dtset%ntypat,dtset%ntypalch,pspheads)
298 
299  psps%npspalch = dtset%npspalch
300  psps%ntypat   = dtset%ntypat
301  psps%ntypalch = dtset%ntypalch
302  psps%ntyppure = dtset%ntyppure
303 
304 !Set the flag for reciprocal space or real space calculations
305  psps%vlspl_recipSpace = (dtset%icoulomb /= 1)
306  psps%positron = dtset%positron
307  psps%useylm   = dtset%useylm
308  psps%usewvl   = dtset%usewvl
309 
310 ! Define treatment of the model core density for NC pseudos.
311  psps%nc_xccc_gspace = dtset%nc_xccc_gspace
312 
313  if (idtset > 1) then
314    ABI_SFREE(psps%algalch)
315    ABI_SFREE(psps%mixalch)
316  end if
317 
318  ABI_MALLOC(psps%algalch,(psps%ntypalch))
319  ABI_MALLOC(psps%mixalch,(psps%npspalch,psps%ntypalch))
320  psps%algalch(1:psps%ntypalch)=dtset%algalch(1:psps%ntypalch)
321 !This value will be overwritten elsewhere in case there are different images ...
322  psps%mixalch(1:psps%npspalch,1:psps%ntypalch)=dtset%mixalch_orig(1:psps%npspalch,1:psps%ntypalch,1)
323 
324 !Set mpspso and psps%pspso
325 !Warning: mpspso might be different for each dataset.
326 !         mpspso not relevant in case of PAW.
327  psps%mpspso=1
328  do ipsp=1,dtset%npsp
329    if(dtset%nspinor==1)then
330      psps%pspso(ipsp)=0
331      ! This is needed to treate SOC perturbatively in sigma.
332      !if (dtset%optdriver == RUNL_SIGMA .and. dtset%so_psp(ipsp) /= 0) then
333      !  ABI_WARNING("Setting pspso to 2 although nspinor == 1")
334      !  psps%pspso(ipsp) = 2
335      !end if
336 
337      ! Ideally the following line should not exist, but at present, the space has to be booked
338      if(pspheads(ipsp)%pspso/=0)psps%mpspso=2
339    else if (psps%usepaw==0) then
340      if(dtset%so_psp(ipsp)/=1)then
341        psps%pspso(ipsp)=dtset%so_psp(ipsp)
342      else
343        psps%pspso(ipsp)=pspheads(ipsp)%pspso
344      end if
345      if(psps%pspso(ipsp)/=0)psps%mpspso=2
346      if(pspheads(ipsp)%pspso/=0)psps%mpspso=2
347    else
348      psps%pspso(ipsp)=1+dtset%pawspnorb
349    end if
350  end do
351 
352 !Set mpssoang, lmnmax, lnmax
353  if(psps%mpspso==1)then
354    psps%mpssoang=psps%mpsang
355    psps%lmnmax  =lmnmax
356    psps%lnmax   =lnmax
357  else
358    psps%mpssoang=2*psps%mpsang-1
359    psps%lmnmax=lmnmaxso
360    psps%lnmax=lnmaxso
361  end if
362 
363 !T. Rangel: for wvl + paw do not change psps%lmnmax
364  if (psps%useylm==0 .and. psps%usepaw/=1 ) then
365    psps%lmnmax=psps%lnmax
366  end if
367 
368 !Set dimekb
369  if (psps%usepaw==0) then
370    psps%dimekb=psps%lnmax
371  else
372    psps%dimekb=psps%lmnmax*(psps%lmnmax+1)/2
373  end if
374 
375 !The following arrays are often not deallocated before the end of the dtset loop
376 !and might keep their content from one dataset to the other, if the conditions are fulfilled
377  changed = 0
378 
379  if(dimekb_old/=psps%dimekb .or. ntypat_old/=dtset%ntypat .or. usepaw_old/=psps%usepaw) then
380    changed = changed + 1
381    if(idtset/=1) then
382      ABI_SFREE(psps%ekb)
383    end if
384    ABI_MALLOC(psps%ekb,(psps%dimekb,dtset%ntypat*(1-psps%usepaw)))
385    psps%ekb = zero
386    dimekb_old=psps%dimekb
387  end if
388 
389  if(lmnmax_old/=psps%lmnmax .or. ntypat_old/=dtset%ntypat)then
390    changed = changed + 1
391    if(idtset/=1) then
392      ABI_SFREE(psps%indlmn)
393    end if
394    ABI_MALLOC(psps%indlmn,(6,psps%lmnmax,dtset%ntypat))
395    psps%indlmn = zero
396    lmnmax_old=psps%lmnmax
397  end if
398 
399  if(mqgridff_old/=psps%mqgrid_ff .or. lnmax_old/=psps%lnmax .or. ntypat_old/=dtset%ntypat)then
400    changed = changed + 1
401    if(idtset/=1) then
402      ABI_SFREE(psps%ffspl)
403      ABI_SFREE(psps%qgrid_ff)
404    end if
405    ABI_MALLOC(psps%ffspl,(psps%mqgrid_ff,2,psps%lnmax,dtset%ntypat))
406    ABI_MALLOC(psps%qgrid_ff,(psps%mqgrid_ff))
407    psps%ffspl = zero
408    psps%qgrid_ff = zero
409    mqgridff_old=psps%mqgrid_ff
410    lnmax_old=psps%lnmax
411  end if
412 
413  if(mqgridvl_old/=psps%mqgrid_vl .or. ntypat_old/=dtset%ntypat)then
414    changed = changed + 1
415    if(idtset/=1) then
416      ABI_SFREE(psps%qgrid_vl)
417      ABI_SFREE(psps%vlspl)
418      if (allocated(psps%nctab)) then
419        do ii=1,size(psps%nctab)
420          call nctab_free(psps%nctab(ii))
421        end do
422        ABI_FREE(psps%nctab)
423      end if
424    end if
425    if (idtset/=1 .and. .not.psps%vlspl_recipSpace) then
426      ABI_SFREE(psps%dvlspl)
427    end if
428 
429    ABI_MALLOC(psps%qgrid_vl,(psps%mqgrid_vl))
430    ABI_MALLOC(psps%vlspl,(psps%mqgrid_vl,2,dtset%ntypat))
431    psps%qgrid_vl = zero
432    psps%vlspl = zero
433 
434    if (psps%usepaw == 0) then
435      ! If you change usepaw in the input, you will get what you deserve!
436      ABI_MALLOC(psps%nctab, (dtset%ntypat))
437      do itypat=1,dtset%ntypat
438        call nctab_init(psps%nctab(itypat), psps%mqgrid_vl, .False., .False.)
439      end do
440    end if
441 
442    if (.not.psps%vlspl_recipSpace) then
443      ABI_MALLOC(psps%dvlspl,(psps%mqgrid_vl,2,dtset%ntypat))
444      psps%dvlspl = zero
445    end if
446    mqgridvl_old=psps%mqgrid_vl
447  end if
448 
449  if(ntypat_old/=dtset%ntypat.or. usepaw_old/=psps%usepaw)then
450    changed = changed + 1
451    if(idtset/=1) then
452      ABI_SFREE(psps%xccc1d)
453    end if
454    ABI_MALLOC(psps%xccc1d,(psps%n1xccc*(1-psps%usepaw),6,dtset%ntypat))
455    psps%xccc1d = zero
456    usepaw_old=psps%usepaw
457  end if
458 
459  if(ntypat_old/=dtset%ntypat)then
460    changed = changed + 1
461    if(idtset/=1) then
462      ABI_SFREE(psps%xcccrc)
463      ABI_SFREE(psps%ziontypat)
464      ABI_SFREE(psps%znucltypat)
465    end if
466    ABI_MALLOC(psps%xcccrc,(dtset%ntypat))
467    ABI_MALLOC(psps%znucltypat,(dtset%ntypat))
468    ABI_MALLOC(psps%ziontypat,(dtset%ntypat))
469    psps%xcccrc = zero
470    psps%znucltypat = zero
471    psps%ziontypat = zero
472    ntypat_old=dtset%ntypat
473  end if
474 
475  psps%ziontypat(:)=dtset%ziontypat(:)
476 
477 end subroutine psps_init_from_dtset

m_psps/psps_init_global [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_init_global

FUNCTION

 Allocate and initialise all part of psps structure that are independent of a given dataset.

INPUTS

 npsp=the number of read pseudo files.
 pspheads(npsp)=<type pspheader_type>all the important information from the
   pseudopotential file header, as well as the psp file name

OUTPUT

SIDE EFFECTS

 psps=<type pseudopotential_type>the pseudopotentials description

SOURCE

176 subroutine psps_init_global(psps, mtypalch, npsp, pspheads)
177 
178 !Arguments ------------------------------------
179 !scalars
180  class(pseudopotential_type),intent(inout) :: psps
181  integer,intent(in) :: mtypalch,npsp
182 !arrays
183  type(pspheader_type),intent(in) :: pspheads(npsp)
184 
185 !Local variables-------------------------------
186 !scalars
187  integer :: ii, mpsang, n1xccc
188 
189 ! *************************************************************************
190 
191 !Allocation of some arrays independent of the dataset
192  ABI_MALLOC(psps%filpsp,(npsp))
193  ABI_MALLOC(psps%pspcod,(npsp))
194  ABI_MALLOC(psps%pspdat,(npsp))
195  ABI_MALLOC(psps%pspso,(npsp))
196  ABI_MALLOC(psps%pspxc,(npsp))
197  ABI_MALLOC(psps%title,(npsp))
198  ABI_MALLOC(psps%zionpsp,(npsp))
199  ABI_MALLOC(psps%znuclpsp,(npsp))
200  call psp2params_init(psps%gth_params, npsp)
201 
202  psps%filpsp(1:npsp)=pspheads(1:npsp)%filpsp
203  psps%pspcod(1:npsp)=pspheads(1:npsp)%pspcod
204  psps%pspdat(1:npsp)=pspheads(1:npsp)%pspdat
205  psps%pspso(1:npsp)=pspheads(1:npsp)%pspso
206  psps%pspxc(1:npsp)=pspheads(1:npsp)%pspxc
207  psps%title(1:npsp)=pspheads(1:npsp)%title
208  psps%zionpsp(1:npsp)=pspheads(1:npsp)%zionpsp
209  psps%znuclpsp(1:npsp)=pspheads(1:npsp)%znuclpsp
210 
211  ! Transfer md5 checksum
212  ABI_MALLOC(psps%md5_pseudos, (npsp))
213  psps%md5_pseudos = pspheads(1:npsp)%md5_checksum
214 !Set values independant from dtset
215  psps%npsp   = npsp
216 !Note that mpsang is the max of 1+lmax, with minimal value 1 (even for local psps, at present)
217  mpsang=1
218  n1xccc=pspheads(1)%xccc
219  do ii=1,psps%npsp
220    mpsang=max(pspheads(ii)%lmax+1,mpsang)
221    n1xccc=max(pspheads(ii)%xccc,n1xccc)
222  end do
223  psps%mpsang = mpsang
224  psps%n1xccc = n1xccc
225 ! Determine here whether the calculation is PAW
226 ! If paw, all pspcod necessarily are 7 or 17 (see iofn2)
227  psps%usepaw  =0
228  if (pspheads(1)%pspcod==7.or.pspheads(1)%pspcod==17) psps%usepaw=1
229  psps%mtypalch = mtypalch
230 
231 end subroutine psps_init_global

m_psps/psps_ncread [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_ncread

FUNCTION

  Read the most important arrays defined in the pseudopotential_type
  in NETCDF file format.
  This function should be called by master node only.

INPUTS

OUTPUT

SOURCE

1113 subroutine psps_ncread(psps, ncid)
1114 
1115 !Arguments ------------------------------------
1116  type(pseudopotential_type),intent(inout) :: psps
1117  integer,intent(in) :: ncid
1118 
1119 !Local variables-------------------------------
1120 !scalars
1121  integer :: ipsp,itypat
1122  integer :: ncerr
1123  integer :: with_xccc
1124 
1125 ! *********************************************************************
1126 
1127  ! Note: Some dimensions and variables are written conditionally,
1128  !       so try to read those but ignore errors
1129 
1130  call psps_free(psps)
1131 
1132  psps%dimekb         = zero
1133  psps%lmnmax         = zero
1134  psps%lnmax          = zero
1135  psps%mproj          = zero
1136  psps%mpsang         = zero
1137  psps%mpspso         = zero
1138  psps%mpssoang       = zero
1139  psps%mqgrid_ff      = zero
1140  psps%mqgrid_vl      = zero
1141  psps%mtypalch       = zero
1142  psps%npsp           = zero
1143  psps%npspalch       = zero
1144  psps%ntypat         = zero
1145  psps%ntypalch       = zero
1146  psps%ntyppure       = zero
1147  psps%n1xccc         = zero
1148  psps%optnlxccc      = zero
1149  psps%positron       = zero
1150  psps%usepaw         = zero
1151  psps%usewvl         = zero
1152  psps%useylm         = zero
1153  psps%nc_xccc_gspace = zero
1154  psps%vlspl_recipSpace = .false.
1155 
1156 #ifdef HAVE_NETCDF
1157 
1158  ! Read dimensions
1159  NCF_CHECK(nctk_get_dim(ncid, "ntypat", psps%ntypat))
1160  NCF_CHECK(nctk_get_dim(ncid, "npsp", psps%npsp))
1161  NCF_CHECK(nctk_get_dim(ncid, "lnmax", psps%lnmax))
1162  NCF_CHECK(nctk_get_dim(ncid, "lmnmax", psps%lmnmax))
1163  NCF_CHECK(nctk_get_dim(ncid, "dimekb", psps%dimekb))
1164  NCF_CHECK(nctk_get_dim(ncid, "mqgrid_vl", psps%mqgrid_vl))
1165  NCF_CHECK(nctk_get_dim(ncid, "mqgrid_ff", psps%mqgrid_ff))
1166  NCF_CHECK(nctk_get_dim(ncid, "n1xccc", psps%n1xccc))
1167  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "usepaw"), psps%usepaw))
1168  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "useylm"), psps%useylm))
1169  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "with_xccc"), with_xccc))
1170 
1171  if (psps%usepaw > 0) with_xccc = 0
1172  if (with_xccc == 0) psps%n1xccc = 0
1173 
1174  ! Allocate arrays
1175  call psp2params_init(psps%gth_params, psps%npsp)
1176  ABI_MALLOC(psps%filpsp,(psps%npsp))
1177  ABI_MALLOC(psps%title,(psps%npsp))
1178  ABI_MALLOC(psps%md5_pseudos,(psps%npsp))
1179 
1180  ABI_MALLOC(psps%pspcod,(psps%npsp))
1181  ABI_MALLOC(psps%pspdat,(psps%npsp))
1182  ABI_MALLOC(psps%pspxc,(psps%npsp))
1183  ABI_MALLOC(psps%pspso,(psps%npsp))
1184 
1185  psps%pspcod = zero
1186  psps%pspdat = zero
1187  psps%pspxc = zero
1188  psps%pspso = zero
1189 
1190  ! GA: zionpsp and znuclpsp dont get written. We assume they are the same
1191  ! as ziontypat and znucltypat
1192  ABI_MALLOC(psps%zionpsp,(psps%npsp))
1193  ABI_MALLOC(psps%znuclpsp,(psps%npsp))
1194  ABI_MALLOC(psps%ziontypat,(psps%ntypat))
1195  ABI_MALLOC(psps%znucltypat,(psps%ntypat))
1196  ABI_MALLOC(psps%xcccrc,(psps%ntypat))
1197  ABI_MALLOC(psps%qgrid_vl,(psps%mqgrid_vl))
1198  ABI_MALLOC(psps%qgrid_ff,(psps%mqgrid_ff))
1199  ABI_MALLOC(psps%indlmn,(6,psps%lmnmax,psps%ntypat))
1200  ABI_MALLOC(psps%vlspl,(psps%mqgrid_vl,2,psps%ntypat))
1201  ABI_MALLOC(psps%ffspl,(psps%mqgrid_ff,2,psps%lmnmax,psps%ntypat))
1202  ABI_MALLOC(psps%ekb,(psps%dimekb,psps%ntypat * (1 - psps%usepaw)))
1203  ABI_MALLOC(psps%xccc1d,(psps%n1xccc,6,psps%ntypat))
1204  ABI_MALLOC(psps%nctab,(psps%ntypat))
1205  if (psps%usepaw == 0) then
1206    do itypat=1,psps%ntypat
1207      psps%nctab(itypat)%mqgrid_vl  = psps%mqgrid_vl
1208      psps%nctab(itypat)%dncdq0     = zero
1209      psps%nctab(itypat)%d2ncdq0    = zero
1210      psps%nctab(itypat)%dnvdq0     = zero
1211      psps%nctab(itypat)%num_tphi   = zero
1212      psps%nctab(itypat)%has_jtot   = .False.
1213 
1214      psps%nctab(itypat)%has_tvale  = .False.
1215      psps%nctab(itypat)%has_tcore  = .False.
1216      ! GA: Do we even need those?
1217      ABI_MALLOC(psps%nctab(itypat)%tvalespl,(psps%mqgrid_vl,2))
1218      ABI_MALLOC(psps%nctab(itypat)%tcorespl,(psps%mqgrid_vl,2))
1219      psps%nctab(itypat)%tvalespl = zero
1220      psps%nctab(itypat)%tcorespl = zero
1221    end do
1222  end if
1223 
1224  psps%ekb = zero
1225  psps%indlmn = zero
1226  psps%xccc1d = zero
1227  psps%xcccrc = zero
1228  psps%vlspl = zero
1229  psps%ffspl = zero
1230  psps%qgrid_vl = zero
1231  psps%qgrid_ff = zero
1232 
1233  psps%zionpsp = zero
1234  psps%znuclpsp = zero
1235  psps%ziontypat = zero
1236  psps%znucltypat = zero
1237 
1238  ! Read variables
1239  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "ziontypat"), psps%ziontypat))
1240  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "znucltypat"), psps%znucltypat))
1241  ! Not dealing with alchemical at the moment.
1242 
1243 
1244  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "spinorbit"), psps%pspso))
1245 
1246  do ipsp=1,psps%npsp
1247    ncerr = nf90_get_var(ncid, nctk_idname(ncid, "filpsp"), psps%filpsp(ipsp), start=[1,ipsp])
1248    ncerr = nf90_get_var(ncid, nctk_idname(ncid, "md5_pseudos"), psps%md5_pseudos(ipsp), start=[1,ipsp])
1249    psps%title(ipsp) = ''
1250  end do
1251  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "qgrid_vl"), psps%qgrid_vl))
1252  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "qgrid_ff"), psps%qgrid_ff))
1253  NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "indlmn"), psps%indlmn))
1254  ncerr = nf90_get_var(ncid, nctk_idname(ncid, "vlspl"), psps%vlspl)
1255  ncerr = nf90_get_var(ncid, nctk_idname(ncid, "ffspl"), psps%ffspl)
1256 
1257  if (psps%usepaw == 0) then
1258    ncerr = nf90_get_var(ncid, nctk_idname(ncid, "ekb"), psps%ekb)
1259 
1260    if (with_xccc > 0) then
1261      NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "xcccrc"), psps%xcccrc))
1262      NCF_CHECK(nf90_get_var(ncid, nctk_idname(ncid, "xccc1d"), psps%xccc1d))
1263    end if
1264 
1265    ! GA: Why bother reading it?
1266    do itypat=1,psps%ntypat
1267      ncerr = nf90_get_var(ncid, nctk_idname(ncid, "nc_tvalespl"), psps%nctab(itypat)%tvalespl, start=[1,1,itypat])
1268      ncerr = nf90_get_var(ncid, nctk_idname(ncid, "nc_tcorespl"), psps%nctab(itypat)%tcorespl, start=[1,1,itypat])
1269    end do
1270 
1271  end if
1272 
1273 #else
1274  ABI_WARNING("netcdf support not activated. psps file cannot be read!")
1275 #endif
1276 
1277 end subroutine psps_ncread

m_psps/psps_ncwrite [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_ncwrite

FUNCTION

  Output the most important arrays defined in the pseudopotential_type
  in NETCDF file format for futher post-processing.
  This function should be called by master node only.

INPUTS

  ncid=NC file handle.

SOURCE

 912 subroutine psps_ncwrite(psps, ncid)
 913 
 914 !Arguments ------------------------------------
 915 !scalars
 916  integer,intent(in) :: ncid
 917  type(pseudopotential_type),intent(in) :: psps
 918 
 919 !Local variables-------------------------------
 920 !scalars
 921  integer :: ipsp,itypat,ncerr
 922  integer :: with_xccc, n1xccc, with_alch
 923 !arrays
 924  real(dp), allocatable :: dummy3(:,:,:)
 925  !real(dp), allocatable :: dummy1(:)
 926 
 927 ! *************************************************************************
 928 
 929 #ifdef HAVE_NETCDF
 930 
 931  with_alch = 0  ! Alchemical IO not supported at the moment.
 932  !psps%mtypalch = zero
 933 
 934  ! GA: Note that lnmax is not used in the DDB text format,
 935  !     so lnmax and lmnmax may be inconsistent in the netcdf file.
 936  !NCF_CHECK(nctk_set_defmode(ncid))
 937 
 938  with_xccc = 0
 939  if (psps%n1xccc > 0) with_xccc = 1
 940  n1xccc = max(1, psps%n1xccc)
 941 
 942  if (.not. allocated(psps%xcccrc)) with_xccc = 0
 943  if (.not. allocated(psps%xccc1d)) with_xccc = 0
 944  if (psps%usepaw /= 0) with_xccc = 0
 945 
 946  ! Define dimensions
 947  ncerr = nctk_def_dims(ncid, [ &
 948      nctkdim_t("fnlen", fnlen + 1), &
 949      nctkdim_t("md5_slen", md5_slen + 1), &
 950      nctkdim_t("ntypat", psps%ntypat), &
 951      nctkdim_t("npsp", psps%npsp), &
 952      nctkdim_t("lnmax", psps%lnmax), &
 953      nctkdim_t("lmnmax", psps%lmnmax), &
 954      nctkdim_t("dimekb", psps%dimekb), &
 955      nctkdim_t("mqgrid_vl", psps%mqgrid_vl), &
 956      nctkdim_t("mqgrid_ff", psps%mqgrid_ff), &
 957      nctkdim_t("n1xccc", n1xccc) &
 958  ])
 959  NCF_CHECK(ncerr)
 960 
 961  ! Define variables
 962  ncerr = nctk_def_iscalars(ncid, [character(len=nctk_slen) :: &
 963                                   "usepaw", "useylm", "with_xccc", "with_alch"])
 964  NCF_CHECK(ncerr)
 965 
 966  ! Arrays
 967  ncerr = nctk_def_arrays(ncid, [&
 968    nctkarr_t("ziontypat", "dp", "ntypat"), &
 969    nctkarr_t("znucltypat", "dp", "ntypat"), &
 970    nctkarr_t("spinorbit", "int", "npsp"), &
 971    nctkarr_t("qgrid_vl", "dp", "mqgrid_vl"), &
 972    nctkarr_t("qgrid_ff", "dp", "mqgrid_ff"), &
 973    nctkarr_t("vlspl", "dp", "mqgrid_vl, two, ntypat"), &
 974    nctkarr_t("indlmn", "int", "six, lmnmax, ntypat"), &
 975    nctkarr_t("ffspl", "dp", "mqgrid_ff, two, lnmax, ntypat"), &
 976    nctkarr_t("filpsp", "char", "fnlen, npsp"), &
 977    nctkarr_t("md5_pseudos", "char", "md5_slen, npsp") &
 978  ])
 979  NCF_CHECK(ncerr)
 980 
 981  if (psps%usepaw == 0) then
 982    NCF_CHECK(nctk_def_arrays(ncid, nctkarr_t("ekb", "dp", "dimekb, ntypat")))
 983    !if (with_xccc > 0) then
 984    NCF_CHECK(nctk_def_arrays(ncid, nctkarr_t("xccc1d", "dp", "n1xccc, six, ntypat")))
 985    NCF_CHECK(nctk_def_arrays(ncid, nctkarr_t("xcccrc", "dp", "ntypat")))
 986 
 987    ncerr = nctk_def_arrays(ncid, [&
 988      nctkarr_t("nc_tvalespl", "dp", "mqgrid_vl, two, ntypat"), &
 989      nctkarr_t("nc_tcorespl", "dp", "mqgrid_vl, two, ntypat")  &
 990    ])
 991    NCF_CHECK(ncerr)
 992  end if
 993 
 994  ! Write data
 995  NCF_CHECK(nf90_put_var(ncid, vid("ziontypat"), psps%ziontypat))
 996  NCF_CHECK(nf90_put_var(ncid, vid("znucltypat"), psps%znucltypat))
 997  ! Note that znuclpsp and ziopsp are not read, since we set with_alch=0
 998 
 999  ncerr = nctk_write_iscalars(ncid, [character(len=nctk_slen) :: &
1000                              "usepaw", "useylm", "with_xccc", "with_alch"], &
1001                              [psps%usepaw, psps%useylm, with_xccc, with_alch])
1002  NCF_CHECK(ncerr)
1003 
1004  if (allocated(psps%pspso)) then
1005     NCF_CHECK(nf90_put_var(ncid, vid("spinorbit"), psps%pspso))
1006  end if
1007  do ipsp=1,psps%npsp
1008    NCF_CHECK(nf90_put_var(ncid, vid("filpsp"), trim(psps%filpsp(ipsp)), start=[1, ipsp]))
1009    NCF_CHECK(nf90_put_var(ncid, vid("md5_pseudos"), trim(psps%md5_pseudos(ipsp)), start=[1, ipsp]))
1010  end do
1011  if (allocated(psps%qgrid_vl)) then
1012    NCF_CHECK(nf90_put_var(ncid, vid("qgrid_vl"), psps%qgrid_vl))
1013  end if
1014  if (allocated(psps%qgrid_ff)) then
1015    NCF_CHECK(nf90_put_var(ncid, vid("qgrid_ff"), psps%qgrid_ff))
1016  end if
1017  if (allocated(psps%indlmn)) then
1018    NCF_CHECK(nf90_put_var(ncid, vid("indlmn"), psps%indlmn))
1019  end if
1020 
1021  ! Local part in q-space and second derivative
1022  if (allocated(psps%vlspl)) then
1023    NCF_CHECK(nf90_put_var(ncid, vid("vlspl"), psps%vlspl))
1024  end if
1025 
1026  ! Form factors for each type of atom
1027  ! for each type and each (l,n) channel, ffnl(q) and second derivative
1028  if (allocated(psps%ffspl)) then
1029    NCF_CHECK(nf90_put_var(ncid, vid("ffspl"), psps%ffspl))
1030  end if
1031 
1032  if (with_xccc > 0) then
1033 
1034  ! Pseudo-core charge for each type of atom, on the real-space radial
1035    NCF_CHECK(nf90_put_var(ncid, vid("xcccrc"), psps%xcccrc))
1036    NCF_CHECK(nf90_put_var(ncid, vid("xccc1d"), psps%xccc1d))
1037 
1038  !else
1039 
1040  !  ABI_MALLOC(dummy1, (psps%ntypat))
1041  !  dummy1 = zero
1042  !  NCF_CHECK(nf90_put_var(ncid, vid("xcccrc"), dummy1))
1043  !  ABI_FREE(dummy1)
1044 
1045  !  ABI_MALLOC(dummy3, (n1xccc, 6, psps%ntypat))
1046  !  dummy3 = zero
1047  !  NCF_CHECK(nf90_put_var(ncid, vid("xccc1d"), dummy3))
1048  !  ABI_FREE(dummy3)
1049 
1050  end if
1051 
1052  ! NC-only: add tcore_spl and tvalespl in q-space
1053  if (psps%usepaw == 0) then
1054    if (allocated(psps%ekb)) then
1055      NCF_CHECK(nf90_put_var(ncid, vid("ekb"), psps%ekb))
1056    end if
1057    do itypat=1,psps%ntypat
1058 
1059      ! TODO Could write variables has_tvale and has_tcore
1060      if (psps%nctab(itypat)%has_tvale) then
1061        ncerr = nf90_put_var(ncid, vid("nc_tvalespl"), psps%nctab(itypat)%tvalespl, start=[1,1,itypat])
1062        NCF_CHECK(ncerr)
1063      else
1064        ABI_MALLOC(dummy3, (psps%mqgrid_vl, 2, psps%ntypat))
1065        dummy3 = zero
1066        ncerr = nf90_put_var(ncid, vid("nc_tvalespl"), dummy3)
1067        NCF_CHECK(ncerr)
1068        ABI_FREE(dummy3)
1069      end if
1070      if (psps%nctab(itypat)%has_tcore) then
1071        ncerr = nf90_put_var(ncid, vid("nc_tcorespl"), psps%nctab(itypat)%tcorespl, start=[1,1,itypat])
1072        NCF_CHECK(ncerr)
1073      else
1074        ABI_MALLOC(dummy3, (psps%mqgrid_vl, 2, psps%ntypat))
1075        dummy3 = zero
1076        ncerr = nf90_put_var(ncid, vid("nc_tcorespl"), dummy3)
1077        NCF_CHECK(ncerr)
1078        ABI_FREE(dummy3)
1079      end if
1080    end do
1081  end if
1082 
1083 #else
1084  ABI_WARNING("netcdf support not activated. psps file cannot be created!")
1085 #endif
1086 
1087 contains
1088  integer function vid(vname)
1089    character(len=*),intent(in) :: vname
1090    vid = nctk_idname(ncid, vname)
1091  end function vid
1092 
1093 end subroutine psps_ncwrite

m_psps/psps_ncwrite_path [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_ncwrite_path

FUNCTION

  Create a new NETCDF file,
  and output the most important arrays defined in the pseudopotential_type
  for futher post-processing.
  This function should be called by master node only.

INPUTS

   path=File name.

SOURCE

871 subroutine psps_ncwrite_path(psps, path)
872 
873 !Arguments ------------------------------------
874  class(pseudopotential_type),intent(in) :: psps
875  character(len=*),intent(in) :: path
876 
877 !Local variables-------------------------------
878  integer :: ncid
879 
880 ! *************************************************************************
881 
882 #ifdef HAVE_NETCDF
883  NCF_CHECK(nctk_open_create(ncid, path, xmpi_comm_self))
884 
885  call psps_ncwrite(psps, ncid)
886 
887  NCF_CHECK(nf90_close(ncid))
888 
889 #else
890  ABI_WARNING("netcdf support not activated. psps file cannot be created!")
891 #endif
892 
893 end subroutine psps_ncwrite_path

m_psps/psps_print [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

 psps_print

FUNCTION

  Print the content of a pseudopotential_type derived type

INPUTS

  psps=<type pseudopotential_type>=Info on the pseudopotentials.
  unit(optional)=unit number for output
  prtvol(optional)=verbosity level
  mode_paral(optional): either "COLL" or "PERS"

OUTPUT

  Only writing

SOURCE

668 subroutine psps_print(psps, unit, prtvol, mode_paral)
669 
670 !Arguments ------------------------------------
671 !scalars
672  class(pseudopotential_type),intent(in) :: psps
673  integer,intent(in),optional :: prtvol,unit
674  character(len=4),intent(in),optional :: mode_paral
675 
676 !Local variables-------------------------------
677 !scalars
678  integer :: ierr,ips,ipsp_alch,ityp_alch,itypat,unt,my_prtvol
679  character(len=4) :: mode
680  character(len=500) :: msg
681 !arrays
682  integer :: cond_values(4)
683  character(len=9) :: cond_string(4)
684 
685 ! *************************************************************************
686 
687  ! Provide defaults
688  my_prtvol=0; if (present(prtvol)) my_prtvol=prtvol
689  unt=std_out; if (present(unit)) unt=unit
690  mode='COLL'; if (present(mode_paral)) mode=mode_paral
691  ierr=0; cond_string(1:4)=' '; cond_values(:)=0
692 
693  ! General info including spin-orbit
694  call wrtout(unt,' ==== Info on pseudopotentials ==== ', mode)
695 
696  SELECT CASE (psps%usepaw)
697  CASE (0)
698    call wrtout(unt,'  Norm-conserving pseudopotentials ', mode)
699    !call wrtout(unt, sjoin('  Max number of Kleinman-Bylander energies ', itoa(psps%dimekb)), mode)
700    !do itypat=1,psps%ntypat
701    ! write(msg,'(a,i4,a,f9.4)')' Type ',itypat,' K-B energies ',(psps%ekb(ikbe,itypat),ikbe=1,psps%dimekb)
702    !end do
703  CASE (1)
704    write(msg,'(a)')
705    call wrtout(unt,'  PAW calculation', mode)
706    !call wrtout(unt,sjoin('  Max number of D_ij coefficients ', itoa(psps%dimekb)), mode)
707  CASE DEFAULT
708    call chkint_eq(0,0,cond_string,cond_values,ierr,'usepaw',psps%usepaw,2,[0,1],unt)
709  END SELECT
710 
711  !SELECT CASE (psps%positron)
712  !CASE (0)
713  !  call wrtout(unt, '  Standard Electron Calculation ', mode)
714  !CASE (1,2)
715  !  write(msg,'(a,i0)')'  Positron Calculation with positron .. ',psps%positron
716  !  call wrtout(unt,msg,mode)
717  !CASE DEFAULT
718  !  call chkint_eq(0,0,cond_string,cond_values,ierr,'positron',psps%positron,3,[0,1,2],unt)
719  !END SELECT
720 
721  write(msg,'(a,i4,2a,i4)')&
722   '  Number of pseudopotentials .. ',psps%npsp,ch10,&
723   '  Number of types of atoms   .. ',psps%ntypat
724  call wrtout(unt,msg,mode)
725 
726  if (psps%usepaw==0) then
727    SELECT CASE (psps%mpspso)
728    CASE (1)
729      call wrtout(unt,'  Scalar calculation (no spin-orbit term) ',mode)
730    CASE (2)
731      write(msg,'(3a,i3)')&
732       '  Calculation with spin-orbit coupling ',ch10,&
733       '  Max number of channels (spin-orbit included) ',psps%mpssoang
734      call wrtout(unt,msg,mode)
735      do itypat=1,psps%ntypat
736        if (psps%pspso(itypat) /= 1) then
737          write(msg,'(a,i4,a,i2,a)')&
738           '  - Atom type ',itypat,' has spin-orbit characteristics (pspso= ',psps%pspso(itypat),")"
739          call wrtout(unt,msg,mode)
740        end if
741      end do
742    CASE DEFAULT
743      call chkint_eq(0,0,cond_string,cond_values,ierr,'mpspso',psps%mpspso,2,[1,2],unt)
744    END SELECT
745  else
746    SELECT CASE (maxval(psps%pspso))
747    CASE (0,1)
748      msg='  Scalar calculation (no spin-orbit term) '
749    CASE (2)
750      msg='  Calculation with spin-orbit coupling '
751    END SELECT
752    call wrtout(unt,msg,mode)
753  end if
754 
755  ! Info on nonlocal part
756  SELECT CASE (psps%useylm)
757  CASE (0)
758    msg = '  Nonlocal part applied using Legendre polynomials '
759  CASE (1)
760    msg = '  Nonlocal part applied using real spherical harmonics '
761  CASE DEFAULT
762    call chkint_eq(0,0,cond_string,cond_values,ierr,'psps%useylm',psps%useylm,2,(/0,1/),unt)
763  END SELECT
764  call wrtout(unt,msg,mode)
765 
766  write(msg,'(a,i3)')'  Max number of non-local projectors over l and type ',psps%mproj
767  call wrtout(unt,msg,mode)
768 
769  write(msg,'(a,i3,2a,i3,2a,i3)')&
770  '  Highest angular momentum +1 ....... ',psps%mpsang,ch10,&
771  '  Max number of (l,n)   components .. ',psps%lnmax, ch10,&
772  '  Max number of (l,m,n) components .. ',psps%lmnmax
773  call wrtout(unt,msg,mode)
774 
775  !FIXME for paw n1xccc==1
776  ! Non-linear Core correction
777  if (psps%n1xccc/=0) then
778    write(msg,'(3a,2(a,i4,a),2a)')ch10,&
779     ' Pseudo-Core Charge Info: ',ch10,&
780     '   Number of radial points for pseudo-core charge .. ',psps%n1xccc,ch10,&
781     '   XC core-correction treatment (optnlxccc) ........ ',psps%optnlxccc,ch10,&
782     '   Radius for pseudo-core charge for each type ..... ',ch10
783    call wrtout(unt,msg,mode)
784    do itypat=1,psps%ntypat
785      write(msg,'(a,i4,a,f7.4)')'  - Atom type ',itypat,' has pseudo-core radius .. ',psps%xcccrc(itypat)
786      call wrtout(unt,msg,mode)
787    end do
788  end if
789 
790  ! Alchemical mixing
791  if (psps%mtypalch/=0) then
792    write(msg,'(3a,3(a,i4,a))')ch10,&
793     ' Calculation with alchemical mixing:',ch10,&
794     '   Number of pure pseudoatoms .... ',psps%ntyppure,ch10,&
795     '   Number of pseudos for mixing .. ',psps%npspalch,ch10,&
796     '   Alchemical pseudoatoms ........ ',psps%ntypalch,ch10
797    call wrtout(unt,msg,mode)
798    do ipsp_alch=1,psps%npspalch
799      do ityp_alch=1,psps%ntypalch
800        write(std_out,*)' mixalch ',psps%mixalch(ipsp_alch,ityp_alch)
801      end do
802    end do
803    do ityp_alch=1,psps%ntypalch
804      write(msg,'(a,i4,a,i4)')' For alchemical atom no. ',ityp_alch,' algalch is .. ',psps%algalch(ityp_alch)
805      call wrtout(unt,msg,mode)
806    end do
807  end if
808 
809  ! Info in Q-grid for spline of form factors
810  write(msg,'(3a,a,i6,a,a,i6)')ch10,&
811   ' Info on the Q-grid used for form factors in spline form: ',ch10,&
812   '   Number of q-points for radial functions ffspl .. ',psps%mqgrid_ff,ch10,&
813   '   Number of q-points for vlspl ................... ',psps%mqgrid_vl
814  call wrtout(unt,msg,mode)
815 
816  if (psps%vlspl_recipSpace) then
817    call wrtout(unt,'   vloc is computed in Reciprocal Space ',mode)
818  else
819    call wrtout(unt,'   vloc is computed in Real Space ',mode)
820  end if
821  if (psps%usepaw == 0) then
822    if (psps%nc_xccc_gspace == 0) call wrtout(unt,'   model core charge treated in real-space', mode)
823    if (psps%nc_xccc_gspace == 1) call wrtout(unt,'   model core charge treated in G-space', mode)
824  end if
825 
826  !TODO additional stuff that might be printed
827  call wrtout(unt, "", mode)
828  do itypat=1,psps%ntypat
829    write(msg,'(a,i0,a,i0)')'  XC functional for type ',itypat,' is ',psps%pspxc(itypat)
830    call wrtout(unt,msg,mode)
831    !write(std_out,*)psps%ziontypat(itypat),psps%znucltypat(itypat)
832    if (psps%usepaw == 0) then
833      call wrtout(unt, sjoin("  Pseudo valence available: ", yesno(psps%nctab(itypat)%has_tvale)), mode)
834    end if
835  end do
836 
837  !integer, pointer :: pspxc(:)
838  ! pspxc(ntypat)
839  ! For each type of psp, the XC functional that was used to generate it, as given by the psp file
840  if (my_prtvol>=3) then
841    do ips=1,psps%npsp
842      write(std_out,*)' Pseudo number   ',ips,' read from ',trim(psps%filpsp(ips))
843      write(std_out,*)' Format or code  ',psps%pspcod(ips)
844      write(std_out,*)' Generation date ',psps%pspdat(ips)
845      write(std_out,*)' Content of first line: ', trim(psps%title(ips))
846    end do
847  end if
848 
849  call wrtout(unt, "", mode)
850 
851 end subroutine psps_print

m_psps/test_xml_xmlpaw_upf [ Functions ]

[ Top ] [ m_psps ] [ Functions ]

NAME

  test_xml_xmlpaw_upf

FUNCTION

  Test if a pseudo potential file is in XML, XML-PAW or in UPF format.

INPUTS

  path=Pseudopotential file

OUTPUT

  usexml=1 if XML file
  xmlpaw=1 if PAW file in XML format
  useupf=1 or 2 if UPF file.

SOURCE

 89 subroutine test_xml_xmlpaw_upf(path, usexml, xmlpaw, useupf)
 90 
 91 !Arguments ------------------------------------
 92 !scalars
 93  character(len=*),intent(in) :: path
 94  integer,intent(out) :: usexml, xmlpaw, useupf
 95 
 96 !Local variables-------------------------------
 97 !scalars
 98  integer :: temp_unit, ii
 99  character(len=500) :: msg,errmsg
100  character(len=70) :: testxml
101 
102 ! *************************************************************************
103 
104 !  Check if the file pseudopotential file is written in XML
105  usexml = 0; xmlpaw = 0; useupf = 0
106 
107  if (open_file(path,msg,newunit=temp_unit,form='formatted',status='old') /= 0) then
108    ABI_ERROR(msg)
109  end if
110  rewind (unit=temp_unit,err=10,iomsg=errmsg)
111 
112  read(temp_unit, "(a)",err=10,iomsg=errmsg) testxml
113  if(testxml(1:5)=='<?xml')then
114    usexml = 1
115    read(temp_unit,*,err=10,iomsg=errmsg) testxml
116    if(testxml(1:4)=='<paw') xmlpaw = 1
117  else
118    usexml = 0
119    if (testxml(1:4) == '<UPF') then
120      ! Make sure this is not UPF version >= 2
121      ! "<UPF version="2.0.1">
122      ii = index(testxml, '"')
123      if (ii /= 0) then
124        useupf = atoi(testxml(ii+1:ii+1))
125        !if (useupf >= 2) then
126        !  ABI_ERROR(sjoin("UPF version >= 2 is not supported by Abinit. Use psp8 or psml format.", ch10, "Pseudo:", path))
127        !end if
128      else
129        ABI_ERROR(sjoin("Cannot find version attributed in UPF file:", path))
130      end if
131 
132    end if
133  end if
134 
135  ! Check if pseudopotential file is a Q-espresso UPF1 file
136  if (useupf == 0) then
137    rewind (unit=temp_unit,err=10,iomsg=errmsg)
138    read(temp_unit,*,err=10,iomsg=errmsg) testxml ! just a string, no relation to xml.
139    if(testxml(1:9)=='<PP_INFO>')then
140      useupf = 1
141    else
142      useupf = 0
143    end if
144  end if
145 
146  close(unit=temp_unit,err=10,iomsg=errmsg)
147 
148  return
149 
150  ! Handle IO error
151 10 continue
152  ABI_ERROR(errmsg)
153 
154 end subroutine test_xml_xmlpaw_upf