TABLE OF CONTENTS


ABINIT/m_pawtab [ Modules ]

[ Top ] [ Modules ]

NAME

  m_pawtab

FUNCTION

  This module contains the definition of the pawtab_type structured datatype,
  as well as related functions and methods.
  pawtab_type variables define TABulated data for PAW (from pseudopotential)

COPYRIGHT

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

NOTES

  FOR DEVELOPPERS: in order to preserve the portability of libPAW library,
  please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt

SOURCE

22 #include "libpaw.h"
23 
24 MODULE m_pawtab
25 
26  USE_DEFS
27  USE_MSG_HANDLING
28  USE_MPI_WRAPPERS
29  USE_MEMORY_PROFILING
30 
31  implicit none
32 
33  private

m_pawtab/pawtab_bcast [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

 pawtab_bcast

FUNCTION

 Communicate pawtab data to all processors

INPUTS

 comm_mpi= communicator used to broadcast data
 [only_from_file]= (optional, default=FALSE)
    If true, only data obtained at the level of the reading
    of the PAW dataset file are broadcasted

SIDE EFFECTS

  pawtab=<type pawtab_type>=a pawtab datastructure

SOURCE

1434 subroutine pawtab_bcast(pawtab,comm_mpi,only_from_file)
1435 
1436 !Arguments ------------------------------------
1437 !scalars
1438  integer,intent(in) :: comm_mpi
1439  logical,optional,intent(in) :: only_from_file
1440  type(pawtab_type),intent(inout) :: pawtab
1441 
1442 !Local variables-------------------------------
1443 !scalars
1444  integer :: ierr,ii,me,nn_dpr,nn_dpr_arr,nn_int,nn_int_arr
1445  integer :: siz_indklmn,siz_indlmn,siz_klmntomn,siz_kmix,siz_lnproju,siz_orbitals
1446  integer :: siz_coredens,siz_coretau,siz_dij0,siz_dltij,siz_dshpfunc,siz_eijkl,siz_eijkl_sr
1447  integer :: siz_euijkl,siz_euij_fll,siz_fk,siz_gammaij,siz_gnorm,siz_fock,siz_kij
1448  integer :: siz_nabla_ij,siz_nabla_im_ij,siz_nablaphi,siz_phi,siz_phiphj,siz_phiphjint
1449  integer :: siz_ph0phiint,siz_qgrid_shp,siz_qijl,siz_rad_for_spline,siz_rhoij0
1450  integer :: siz_shape_alpha,siz_shape_q,siz_shapefunc,siz_shapefncg,siz_sij,siz_tcoredens
1451  integer :: siz_tcoretau,siz_tcorespl,siz_tcoretauspl,siz_tnablaphi,siz_tphi,siz_tphitphj
1452  integer :: siz_tproj,siz_tvalespl,siz_vee,siz_vex,siz_vhtnzc,siz_vhnzc,siz_vminushalf
1453  integer :: siz_zioneff,siz_wvlpaw,siz_wvl_pngau,siz_wvl_parg,siz_wvl_pfac
1454  integer :: siz_wvl_rholoc_rad,siz_wvl_rholoc_d,sz1,sz2
1455  logical :: full_broadcast
1456  character (len=500) :: msg,msg0
1457 !arrays
1458  integer :: nn(4)
1459  integer,allocatable :: list_int(:)
1460  real(dp),allocatable :: list_dpr(:)
1461 
1462 !*************************************************************************
1463 
1464  me=xmpi_comm_rank(comm_mpi)
1465  full_broadcast=.true.;if (present(only_from_file)) full_broadcast=(.not.only_from_file)
1466 
1467  nn_int=0 ; nn_int_arr=0 ; nn_dpr=0 ; nn_dpr_arr=0
1468 
1469 !=========================================================================
1470 !Compute the amount of data to communicate
1471 !=========================================================================
1472 
1473  if (me==0) then
1474    msg=''
1475 
1476 !Integers (read from psp file)
1477 !-------------------------------------------------------------------------
1478 !  basis_size,has_coretau,has_fock,has_kij,has_shapefncg,has_nabla,has_nablaphi,has_tproj
1479 !  has_tvale,has_vhtnzc,has_vhnzc,has_vminushalf,has_wvl,ij_size,l_size,lmn_size,lmn2_size
1480 !  mesh_size,partialwave_mesh_size,core_mesh_size,coretau_mesh_size,vminus_mesh_size
1481 !  tnvale_mesh_size,mqgrid,shape_lambda,shape_type,usetcore,usexcnhat
1482    nn_int=nn_int+28
1483 
1484 !Integers (depending on the parameters of the calculation)
1485 !-------------------------------------------------------------------------
1486 !  ij_proj,lcut_size,lexexch,lmnmix_sz,lpawu,mqgrid_shp,nproju,useexexch,usepawu,usepotzero,
1487 !  option_interaction_pawu,usespnorb
1488    if (full_broadcast) nn_int=nn_int+12
1489 
1490 !Reals (read from psp file)
1491 !-------------------------------------------------------------------------
1492 !  beta,dncdq0,d2ncdq0,dnvdq0,dtaucdq0,ex_cc,exccore,lamb_shielding,rpaw,rshp,rcore,rcoretau,shape_sigma
1493    nn_dpr=nn_dpr+13
1494 
1495 !Reals (depending on the parameters of the calculation)
1496 !-------------------------------------------------------------------------
1497 !  exchmix,f4of2_sla,f6of2_sla,jpawu,upawu
1498    if (full_broadcast) nn_dpr=nn_dpr+5
1499 
1500 !Integers arrays (read from psp file)
1501 !-------------------------------------------------------------------------
1502    siz_indlmn=0 ; siz_orbitals=0
1503    nn_int=nn_int+2
1504    if (allocated(pawtab%indlmn)) then
1505      siz_indlmn=size(pawtab%indlmn)                  !(6,lmn_size)
1506      if (siz_indlmn/=6*pawtab%lmn_size) msg=trim(msg)//' indlmn'
1507      nn_int_arr=nn_int_arr+siz_indlmn
1508    end if
1509    if (allocated(pawtab%orbitals)) then
1510      siz_orbitals=size(pawtab%orbitals)              !(basis_size)
1511      if (siz_orbitals/=pawtab%basis_size) msg=trim(msg)//' orbitals'
1512      nn_int_arr=nn_int_arr+siz_orbitals
1513    end if
1514 
1515 !Integers arrays (depending on the parameters of the calculation)
1516 !-------------------------------------------------------------------------
1517    siz_indklmn=0 ; siz_klmntomn=0 ; siz_kmix=0 ; siz_lnproju=0
1518    if (full_broadcast) then
1519      nn_int=nn_int+4
1520      if (allocated(pawtab%indklmn)) then
1521        siz_indklmn=size(pawtab%indklmn)             !(6,lmn2_size)
1522        if (siz_indklmn/=8*pawtab%lmn2_size) msg=trim(msg)//' indklmn'
1523        nn_int_arr=nn_int_arr+siz_indklmn
1524      end if
1525      if (allocated(pawtab%klmntomn)) then
1526        siz_klmntomn=size(pawtab%klmntomn)           !(4,lmn2_size)
1527        if (siz_klmntomn/=4*pawtab%lmn2_size) msg=trim(msg)//' klmntomn'
1528        nn_int_arr=nn_int_arr+siz_klmntomn
1529      end if
1530      if (allocated(pawtab%kmix)) then
1531        siz_kmix=size(pawtab%kmix)                   !(lmnmix_sz)
1532        if (siz_kmix/=6*pawtab%lmnmix_sz) msg=trim(msg)//' kmix'
1533        nn_int_arr=nn_int_arr+siz_kmix
1534      end if
1535      if (allocated(pawtab%lnproju)) then
1536        siz_lnproju=size(pawtab%lnproju)             !(nproju)
1537        if (siz_lnproju/=pawtab%nproju) msg=trim(msg)//' lnproju'
1538        nn_int_arr=nn_int_arr+siz_lnproju
1539      end if
1540    end if ! full_broadcast
1541 
1542 !Reals arrays (read from psp file)
1543 !-------------------------------------------------------------------------
1544    siz_coredens=0 ; siz_coretau=0    ; siz_dij0=0     ; siz_kij=0        ; siz_fock=0
1545    siz_phi=0      ; siz_rhoij0=0   ; siz_shape_alpha=0
1546    siz_shape_q=0  ; siz_shapefunc=0  ; siz_tcoredens=0; siz_tcoretau=0
1547    siz_tcorespl=0 ; siz_tcoretauspl=0; siz_tphi=0       ; siz_tproj=0
1548    siz_tvalespl=0 ; siz_vhtnzc=0     ; siz_vhnzc=0    ; siz_vminushalf=0
1549    nn_int=nn_int+20
1550 
1551    if (allocated(pawtab%coredens)) then
1552      siz_coredens=size(pawtab%coredens)             !(core_mesh_size)
1553      if (siz_coredens/=pawtab%core_mesh_size) &
1554 &      msg=trim(msg)//' coredens'
1555      nn_dpr=nn_dpr+siz_coredens
1556    end if
1557    if (allocated(pawtab%coretau)) then
1558      siz_coretau=size(pawtab%coretau)             !(coretau_mesh_size)
1559      if (siz_coretau/=pawtab%coretau_mesh_size) &
1560 &      msg=trim(msg)//' coretau'
1561      nn_dpr=nn_dpr+siz_coretau
1562    end if
1563    if (allocated(pawtab%dij0)) then
1564      siz_dij0=size(pawtab%dij0)                     !(lmn2_size)
1565      if (siz_dij0/=pawtab%lmn2_size) msg=trim(msg)//' dij0'
1566      nn_dpr=nn_dpr+siz_dij0
1567    end if
1568    if (allocated(pawtab%ex_cvij)) then
1569      siz_fock=size(pawtab%ex_cvij)                  !(lmn2_size)
1570      if (siz_fock/=pawtab%lmn2_size) msg=trim(msg)//' fock'
1571      nn_dpr=nn_dpr+siz_fock
1572    end if
1573    if (allocated(pawtab%kij)) then
1574      siz_kij=size(pawtab%kij)                       !(lmn2_size)
1575      if (siz_kij/=pawtab%lmn2_size) msg=trim(msg)//' kij'
1576      nn_dpr=nn_dpr+siz_kij
1577    end if
1578    if (allocated(pawtab%phi)) then
1579      siz_phi=size(pawtab%phi)                       !(partialwave_mesh_size, basis_size)
1580      if (siz_phi/=pawtab%partialwave_mesh_size*pawtab%basis_size) msg=trim(msg)//' phi'
1581      nn_dpr=nn_dpr+siz_phi
1582    end if
1583    if (allocated(pawtab%rhoij0)) then
1584      siz_rhoij0=size(pawtab%rhoij0)                 !(lmn2_size)
1585      if (siz_rhoij0/=pawtab%lmn2_size) msg=trim(msg)//' rhoij0'
1586      nn_dpr=nn_dpr+siz_rhoij0
1587    end if
1588    if (allocated(pawtab%shape_alpha)) then
1589      siz_shape_alpha=size(pawtab%shape_alpha)       !(2,l_size)
1590      if (siz_shape_alpha/=pawtab%l_size*2) msg=trim(msg)//' shape_alpha'
1591      nn_dpr=nn_dpr+siz_shape_alpha
1592    end if
1593    if (allocated(pawtab%shape_q)) then
1594      siz_shape_q=size(pawtab%shape_q)               !(2,l_size)
1595      if (siz_shape_q/=pawtab%l_size*2) msg=trim(msg)//' shape_q'
1596      nn_dpr=nn_dpr+siz_shape_q
1597    end if
1598    if (allocated(pawtab%shapefunc)) then
1599      siz_shapefunc=size(pawtab%shapefunc)           !(mesh_size,l_size)
1600      if (siz_shapefunc/=pawtab%mesh_size*pawtab%l_size) msg=trim(msg)//' shapefunc'
1601      nn_dpr=nn_dpr+siz_shapefunc
1602    end if
1603    if (allocated(pawtab%tcoredens)) then
1604      siz_tcoredens=size(pawtab%tcoredens)           !(core_mesh_size,1 or 6)
1605      if (siz_tcoredens/=pawtab%core_mesh_size.and.siz_tcoredens/=6*pawtab%core_mesh_size) &
1606 &      msg=trim(msg)//' tcoredens'
1607      nn_dpr=nn_dpr+siz_tcoredens
1608    end if
1609    if (allocated(pawtab%tcoretau)) then
1610      siz_tcoretau=size(pawtab%tcoretau)             !(coretau_mesh_size,1)
1611      if (siz_tcoretau/=pawtab%coretau_mesh_size) &
1612 &      msg=trim(msg)//' tcoretau'
1613      nn_dpr=nn_dpr+siz_tcoretau
1614    end if
1615    if (allocated(pawtab%tcorespl)) then
1616      siz_tcorespl=size(pawtab%tcorespl)             !(mqgrid,2)
1617      if (siz_tcorespl/=pawtab%mqgrid*2) msg=trim(msg)//' tcorespl'
1618      nn_dpr=nn_dpr+siz_tcorespl
1619    end if
1620    if (allocated(pawtab%tcoretauspl)) then
1621      siz_tcoretauspl=size(pawtab%tcoretauspl)       !(mqgrid,2)
1622      if (siz_tcoretauspl/=pawtab%mqgrid*2.and.siz_tcoretauspl/=0) msg=trim(msg)//' tcoretauspl'
1623      nn_dpr=nn_dpr+siz_tcoretauspl
1624    end if
1625    if (allocated(pawtab%tphi)) then
1626      siz_tphi=size(pawtab%tphi)                     !(partialwave_mesh_size, basis_size)
1627      if (siz_tphi/=pawtab%partialwave_mesh_size*pawtab%basis_size) msg=trim(msg)//' tphi'
1628      nn_dpr=nn_dpr+siz_tphi
1629    end if
1630    if (allocated(pawtab%tproj)) then
1631      siz_tproj=size(pawtab%tproj)                   !(???,basis_size)
1632     if (mod(siz_tproj,pawtab%basis_size)/=0) msg=trim(msg)//' tproj'
1633      nn_dpr=nn_dpr+siz_tproj
1634    end if
1635    if (allocated(pawtab%tvalespl)) then
1636      siz_tvalespl=size(pawtab%tvalespl)             !(mqgrid or mesh_size or tnvale_mesh_size,2)
1637      if (siz_tvalespl/=2*pawtab%mqgrid.and.siz_tvalespl/=2*pawtab%mesh_size.and. &
1638 &        siz_tvalespl/=2*pawtab%tnvale_mesh_size) msg=trim(msg)//' tvalespl'
1639      nn_dpr=nn_dpr+siz_tvalespl
1640    end if
1641    if (allocated(pawtab%vhtnzc)) then
1642      siz_vhtnzc=size(pawtab%vhtnzc)                 !(mesh_size)
1643      if (siz_vhtnzc/=pawtab%mesh_size) msg=trim(msg)//' vhtnzc'
1644      nn_dpr=nn_dpr+siz_vhtnzc
1645    end if
1646    if (allocated(pawtab%vhnzc)) then
1647      siz_vhnzc=size(pawtab%vhnzc)                   !(mesh_size)
1648      if (siz_vhnzc/=pawtab%mesh_size) msg=trim(msg)//' vhnzc'
1649      nn_dpr=nn_dpr+siz_vhnzc
1650    end if
1651    if (allocated(pawtab%vminushalf)) then
1652      siz_vminushalf=size(pawtab%vminushalf)         !(mesh_size)
1653      if (siz_vminushalf/=pawtab%vminus_mesh_size) msg=trim(msg)//' vvminushalf'
1654      nn_dpr=nn_dpr+siz_vminushalf
1655    end if
1656 
1657 !Reals arrays (depending on the parameters of the calculation)
1658 !-------------------------------------------------------------------------
1659    siz_dltij=0    ; siz_dshpfunc=0
1660    siz_eijkl=0    ; siz_eijkl_sr=0 ; siz_euijkl=0    ; siz_euij_fll=0
1661    siz_fk=0       ; siz_gammaij=0  ; siz_gnorm=0
1662    siz_nabla_ij=0 ; siz_nabla_im_ij=0
1663    siz_nablaphi=0 ; siz_phiphj=0   ; siz_phiphjint=0 ; siz_ph0phiint=0
1664    siz_qgrid_shp=0; siz_qijl=0     ; siz_rad_for_spline=0
1665    siz_shapefncg=0; siz_sij=0      ; siz_tnablaphi=0 ; siz_tphitphj=0
1666    siz_vee=0      ; siz_vex=0      ; siz_zioneff=0
1667    if (full_broadcast) then
1668      nn_int=nn_int+25
1669      if (allocated(pawtab%dltij)) then
1670        siz_dltij=size(pawtab%dltij)                   !(lmn2_size)
1671        if (siz_dltij/=pawtab%lmn2_size) msg=trim(msg)//' dltij'
1672        nn_dpr=nn_dpr+siz_dltij
1673      end if
1674      if (allocated(pawtab%dshpfunc)) then
1675        siz_dshpfunc=size(pawtab%dshpfunc)             !(mesh_size,l_size,4)
1676        if (siz_dshpfunc/=pawtab%mesh_size*pawtab%l_size*4) msg=trim(msg)//' dshpfunc'
1677        nn_dpr=nn_dpr+siz_dshpfunc
1678      end if
1679      if (allocated(pawtab%eijkl)) then
1680        siz_eijkl=size(pawtab%eijkl)                   !(lmn2_size,lmn2_size)
1681        if (siz_eijkl/=pawtab%lmn2_size*pawtab%lmn2_size) msg=trim(msg)//' eijkl'
1682        nn_dpr=nn_dpr+siz_eijkl
1683      end if
1684      if (allocated(pawtab%eijkl_sr)) then
1685        siz_eijkl_sr=size(pawtab%eijkl_sr)             !(lmn2_size,lmn2_size)
1686        if (siz_eijkl_sr/=pawtab%lmn2_size*pawtab%lmn2_size) msg=trim(msg)//' eijkl_sr'
1687        nn_dpr=nn_dpr+siz_eijkl_sr
1688      end if
1689      if (allocated(pawtab%euijkl)) then
1690        siz_euijkl=size(pawtab%euijkl)                 !(3,lmn_size,lmn_size,lmn_size,lmn_size)
1691        if (siz_euijkl/=3*pawtab%lmn_size*pawtab%lmn_size*pawtab%lmn_size*pawtab%lmn_size) msg=trim(msg)//' euijkl'
1692        nn_dpr=nn_dpr+siz_euijkl
1693      end if
1694      if (allocated(pawtab%euij_fll)) then
1695        siz_euij_fll=size(pawtab%euij_fll)             !(2,2,lmn_size,lmn_size,lmn_size,lmn_size)
1696        if (siz_euij_fll/=4*pawtab%lmn_size*pawtab%lmn_size*pawtab%lmn_size*pawtab%lmn_size) msg=trim(msg)//' euij_fll'
1697        nn_dpr=nn_dpr+siz_euij_fll
1698      end if
1699      if (allocated(pawtab%fk)) then
1700        siz_fk=size(pawtab%fk)                         !(6,4)
1701        if (siz_fk/=24) msg=trim(msg)//' fk'
1702        nn_dpr=nn_dpr+siz_fk
1703      end if
1704      if (allocated(pawtab%gammaij)) then
1705        siz_gammaij=size(pawtab%gammaij)               !(l_size)
1706        if (siz_gammaij/=pawtab%l_size) msg=trim(msg)//' gammaij'
1707        nn_dpr=nn_dpr+siz_gammaij
1708      end if
1709      if (allocated(pawtab%gnorm)) then
1710        siz_gnorm=size(pawtab%gnorm)                   !(l_size)
1711        if (siz_gnorm/=pawtab%l_size) msg=trim(msg)//' gnorm'
1712        nn_dpr=nn_dpr+siz_gnorm
1713      end if
1714      if (allocated(pawtab%nabla_ij)) then
1715        siz_nabla_ij=size(pawtab%nabla_ij)             !(3,lmn_size,lmn_size)
1716        if (siz_nabla_ij/=pawtab%lmn_size) msg=trim(msg)//' nabla_ij'
1717        nn_dpr=nn_dpr+siz_nabla_ij
1718      end if
1719      if (allocated(pawtab%nabla_im_ij)) then
1720        siz_nabla_im_ij=size(pawtab%nabla_im_ij)       !(3,lmn_size,lmn_size)
1721        if (siz_nabla_im_ij/=pawtab%lmn_size) msg=trim(msg)//' nabla_im_ij'
1722        nn_dpr=nn_dpr+siz_nabla_im_ij
1723      end if
1724      if (allocated(pawtab%nablaphi)) then
1725        siz_phi=size(pawtab%nablaphi)                  !(partialwave_mesh_size, basis_size)
1726        if (siz_nablaphi/=pawtab%partialwave_mesh_size*pawtab%basis_size) msg=trim(msg)//' nablaphi'
1727        nn_dpr=nn_dpr+siz_nablaphi
1728      end if
1729      if (allocated(pawtab%phiphj)) then
1730        siz_phiphj=size(pawtab%phiphj)                 !(mesh_size,ij_size)
1731        if (siz_phiphj/=pawtab%mesh_size*pawtab%ij_size) msg=trim(msg)//' phiphj'
1732        nn_dpr=nn_dpr+siz_phiphj
1733      end if
1734      if (allocated(pawtab%phiphjint)) then
1735        siz_phiphjint=size(pawtab%phiphjint)           !(ij_proj)
1736        if (siz_phiphjint/=pawtab%ij_proj) msg=trim(msg)//' phiphjint'
1737        nn_dpr=nn_dpr+siz_phiphjint
1738      end if
1739      if (allocated(pawtab%ph0phiint)) then
1740        siz_ph0phiint=size(pawtab%ph0phiint)           !(ij_proj)
1741        if (siz_ph0phiint/=pawtab%ij_proj) msg=trim(msg)//' ph0phiint'
1742        nn_dpr=nn_dpr+siz_ph0phiint
1743      end if
1744      if (allocated(pawtab%qgrid_shp)) then
1745        siz_qgrid_shp=size(pawtab%qgrid_shp)           !(mqgrid_shp)
1746        if (siz_qgrid_shp/=pawtab%mqgrid_shp) msg=trim(msg)//' qgrid_shp'
1747        nn_dpr=nn_dpr+siz_qgrid_shp
1748      end if
1749      if (allocated(pawtab%qijl)) then
1750        siz_qijl=size(pawtab%qijl)                     !(l_size**2,lmn2_size)
1751        if (siz_qijl/=pawtab%l_size**2*pawtab%lmn2_size) msg=trim(msg)//' qijl'
1752        nn_dpr=nn_dpr+siz_qijl
1753      end if
1754      if (allocated(pawtab%rad_for_spline)) then
1755        siz_rad_for_spline=size(pawtab%rad_for_spline) !(mesh_size)
1756        if (siz_rad_for_spline/=pawtab%mesh_size) msg=trim(msg)//' rad_for_spline'
1757        nn_dpr=nn_dpr+siz_rad_for_spline
1758      end if
1759      if (allocated(pawtab%shapefncg)) then
1760        siz_shapefncg=size(pawtab%shapefncg)           !(mqgrid_shp,2,l_size)
1761        if (siz_shapefncg/=2*pawtab%mqgrid_shp*pawtab%l_size) msg=trim(msg)//' shapefncg'
1762        nn_dpr=nn_dpr+siz_shapefncg
1763      end if
1764      if (allocated(pawtab%sij)) then
1765        siz_sij=size(pawtab%sij)                       !(lmn2_size)
1766        if (siz_sij/=pawtab%lmn2_size) msg=trim(msg)//' sij'
1767        nn_dpr=nn_dpr+siz_sij
1768      end if
1769      if (allocated(pawtab%tnablaphi)) then
1770        siz_tnablaphi=size(pawtab%tnablaphi)           !(partialwave_mesh_size, basis_size)
1771        if (siz_tnablaphi/=pawtab%partialwave_mesh_size*pawtab%basis_size) msg=trim(msg)//' tnablaphi'
1772        nn_dpr=nn_dpr+siz_tnablaphi
1773      end if
1774      if (allocated(pawtab%tphitphj)) then
1775        siz_tphitphj=size(pawtab%tphitphj)             !(mesh_size,ij_size)
1776        if (siz_tphitphj/=pawtab%mesh_size*pawtab%ij_size) msg=trim(msg)//' tphitphj'
1777        nn_dpr=nn_dpr+siz_tphitphj
1778      end if
1779      if (allocated(pawtab%vee)) then
1780        siz_vee=size(pawtab%vee)                       !(2*lpawu+1,2*lpawu+1,2*lpawu+1,2*lpawu+1)
1781        if (siz_vee/=(2*pawtab%lpawu+1)**4) msg=trim(msg)//' vee'
1782        nn_dpr=nn_dpr+siz_vee
1783      end if
1784      if (allocated(pawtab%vex)) then
1785        siz_vex=size(pawtab%vex)                       !(2*lexexch+1,2*lexexch+1,2*lexexch+1,2*lexexch+1,4)
1786        if (siz_vex/=4*(2*pawtab%lpawu+1)**4) msg=trim(msg)//' vex'
1787        nn_dpr=nn_dpr+siz_vex
1788      end if
1789      if (allocated(pawtab%zioneff)) then
1790        siz_zioneff=size(pawtab%zioneff)               !(ij_proj)
1791        if (siz_zioneff/=pawtab%ij_proj) msg=trim(msg)//' zioneff'
1792        nn_dpr=nn_dpr+siz_zioneff
1793      end if
1794    end if ! full_broadcast
1795 
1796 !Datastructures (read from psp file)
1797 !-------------------------------------------------------------------------
1798    siz_wvl_pngau=0 ; siz_wvl_parg=0 ; siz_wvl_pfac=0
1799    siz_wvl_rholoc_rad=0 ; siz_wvl_rholoc_d=0
1800    siz_wvlpaw=0
1801    nn_int=nn_int+1
1802    if (associated(pawtab%wvl)) then
1803      siz_wvlpaw=1
1804      nn_int=nn_int+3
1805 !    wvl%npspcode_init_guess,wvl%ptotgau
1806      nn_int=nn_int+2
1807      if (allocated(pawtab%wvl%pngau)) then
1808        siz_wvl_pngau=size(pawtab%wvl%pngau)         !(basis_size)
1809        if (siz_wvl_pngau/=pawtab%basis_size) msg=trim(msg)//' wvl_pngau'
1810        nn_int_arr=nn_int_arr+siz_wvl_pngau
1811      end if
1812      if (allocated(pawtab%wvl%parg)) then
1813        siz_wvl_parg=size(pawtab%wvl%parg)          !(2,ptotgau)
1814        if (siz_wvl_parg/=2*pawtab%wvl%ptotgau) msg=trim(msg)//' wvl_parg'
1815        nn_dpr_arr=nn_dpr_arr+siz_wvl_parg
1816      end if
1817      if (allocated(pawtab%wvl%pfac)) then
1818        siz_wvl_pfac=size(pawtab%wvl%pfac )         !(2,ptotgau)
1819        if (siz_wvl_pfac/=2*pawtab%wvl%ptotgau) msg=trim(msg)//' wvl_pfac'
1820        nn_dpr_arr=nn_dpr_arr+siz_wvl_pfac
1821      end if
1822 !    wvl%rholoc%msz
1823      nn_int=nn_int+3
1824      if (pawtab%wvl%rholoc%msz>0) then
1825        if (allocated(pawtab%wvl%rholoc%rad)) then
1826          siz_wvl_rholoc_rad=size(pawtab%wvl%rholoc%rad) !(msz)
1827          if (siz_wvl_rholoc_rad/=pawtab%wvl%rholoc%msz) msg=trim(msg)//' wvl_rholoc_rad'
1828          nn_dpr_arr=nn_dpr_arr+siz_wvl_rholoc_rad
1829        end if
1830        if (allocated(pawtab%wvl%rholoc%d)) then
1831          siz_wvl_rholoc_d=size(pawtab%wvl%rholoc%d)     !(msz,4)
1832          if (siz_wvl_rholoc_d/=4*pawtab%wvl%rholoc%msz) msg=trim(msg)//' wvl_rholoc_d'
1833          nn_dpr_arr=nn_dpr_arr+siz_wvl_rholoc_d
1834        end if
1835      end if
1836    end if
1837 
1838 !Datastructures (depending on the parameters of the calculation)
1839 !-------------------------------------------------------------------------
1840 !  Nothing
1841 
1842 !  Are the sizes OK ?
1843    if (trim(msg)/='') then
1844      write(msg0,'(3a)') &
1845 &     'There is a problem with the size of the following array(s):',ch10,trim(msg)
1846      LIBPAW_BUG(msg0)
1847    end if
1848 
1849  end if ! me=0
1850 
1851 !Broadcast the sizes of buffers
1852 !=========================================================================
1853 
1854  if (me==0) then
1855    nn(1)=nn_int ; nn(2)=nn_int_arr
1856    nn(3)=nn_dpr ; nn(4)=nn_dpr_arr
1857  end if
1858  call xmpi_bcast(nn,0,comm_mpi,ierr)
1859  if (me/=0) then
1860    nn_int=nn(1) ; nn_int_arr=nn(2)
1861    nn_dpr=nn(3) ; nn_dpr_arr=nn(4)
1862  end if
1863 
1864 !Broadcast all the integer: sizes, integer scalars, integer arrays
1865 !=========================================================================
1866 
1867  LIBPAW_ALLOCATE(list_int,(nn_int+nn_int_arr))
1868 
1869 !Fill the buffer of the sender
1870 !-------------------------------------------------------------------------
1871  if (me==0) then
1872    ii=1
1873 
1874 !First the data read from a psp file
1875 !...................................
1876 
1877 !Sizes of arrays (read from psp file)
1878    list_int(ii)=siz_indlmn  ;ii=ii+1
1879    list_int(ii)=siz_orbitals  ;ii=ii+1
1880    list_int(ii)=siz_coredens  ;ii=ii+1
1881    list_int(ii)=siz_coretau  ;ii=ii+1
1882    list_int(ii)=siz_dij0  ;ii=ii+1
1883    list_int(ii)=siz_kij  ;ii=ii+1
1884    list_int(ii)=siz_fock  ;ii=ii+1
1885    list_int(ii)=siz_phi  ;ii=ii+1
1886    list_int(ii)=siz_rhoij0  ;ii=ii+1
1887    list_int(ii)=siz_shape_alpha  ;ii=ii+1
1888    list_int(ii)=siz_shape_q  ;ii=ii+1
1889    list_int(ii)=siz_shapefunc  ;ii=ii+1
1890    list_int(ii)=siz_tcoredens  ;ii=ii+1
1891    list_int(ii)=siz_tcoretau  ;ii=ii+1
1892    list_int(ii)=siz_tcorespl  ;ii=ii+1
1893    list_int(ii)=siz_tcoretauspl  ;ii=ii+1
1894    list_int(ii)=siz_tphi  ;ii=ii+1
1895    list_int(ii)=siz_tproj  ;ii=ii+1
1896    list_int(ii)=siz_tvalespl  ;ii=ii+1
1897    list_int(ii)=siz_vhtnzc  ;ii=ii+1
1898    list_int(ii)=siz_vhnzc  ;ii=ii+1
1899    list_int(ii)=siz_vminushalf  ;ii=ii+1
1900    list_int(ii)=siz_wvlpaw  ;ii=ii+1
1901 !Integers (read from psp file)
1902    list_int(ii)=pawtab%basis_size  ;ii=ii+1
1903    list_int(ii)=pawtab%has_fock  ;ii=ii+1
1904    list_int(ii)=pawtab%has_kij  ;ii=ii+1
1905    list_int(ii)=pawtab%has_shapefncg  ;ii=ii+1
1906    list_int(ii)=pawtab%has_nabla  ;ii=ii+1
1907    list_int(ii)=pawtab%has_nablaphi ; ii=ii+1
1908    list_int(ii)=pawtab%has_tproj  ;ii=ii+1
1909    list_int(ii)=pawtab%has_tvale  ;ii=ii+1
1910    list_int(ii)=pawtab%has_coretau  ;ii=ii+1
1911    list_int(ii)=pawtab%has_vhtnzc  ;ii=ii+1
1912    list_int(ii)=pawtab%has_vhnzc  ;ii=ii+1
1913    list_int(ii)=pawtab%has_vminushalf  ;ii=ii+1
1914    list_int(ii)=pawtab%has_wvl  ;ii=ii+1
1915    list_int(ii)=pawtab%ij_size  ;ii=ii+1
1916    list_int(ii)=pawtab%l_size  ;ii=ii+1
1917    list_int(ii)=pawtab%lmn_size  ;ii=ii+1
1918    list_int(ii)=pawtab%lmn2_size  ;ii=ii+1
1919    list_int(ii)=pawtab%mesh_size  ;ii=ii+1
1920    list_int(ii)=pawtab%partialwave_mesh_size  ;ii=ii+1
1921    list_int(ii)=pawtab%core_mesh_size  ;ii=ii+1
1922    list_int(ii)=pawtab%coretau_mesh_size  ;ii=ii+1
1923    list_int(ii)=pawtab%vminus_mesh_size  ;ii=ii+1
1924    list_int(ii)=pawtab%tnvale_mesh_size  ;ii=ii+1
1925    list_int(ii)=pawtab%mqgrid  ;ii=ii+1
1926    list_int(ii)=pawtab%shape_lambda  ;ii=ii+1
1927    list_int(ii)=pawtab%shape_type  ;ii=ii+1
1928    list_int(ii)=pawtab%usetcore  ;ii=ii+1
1929    list_int(ii)=pawtab%usexcnhat  ;ii=ii+1
1930 !Integer arrays (read from psp file)
1931    if (siz_indlmn>0) then
1932      list_int(ii:ii+siz_indlmn-1)=reshape(pawtab%indlmn,(/siz_indlmn/))
1933      ii=ii+siz_indlmn
1934    end if
1935    if (siz_orbitals>0) then
1936      list_int(ii:ii+siz_orbitals-1)=pawtab%orbitals(1:siz_orbitals)
1937      ii=ii+siz_orbitals
1938    end if
1939 !Integers in datastructures (read from psp file)
1940    if (siz_wvlpaw==1) then
1941      list_int(ii)=siz_wvl_pngau  ;ii=ii+1
1942      list_int(ii)=siz_wvl_parg  ;ii=ii+1
1943      list_int(ii)=siz_wvl_pfac  ;ii=ii+1
1944      list_int(ii)=pawtab%wvl%npspcode_init_guess  ;ii=ii+1
1945      list_int(ii)=pawtab%wvl%ptotgau  ;ii=ii+1
1946      if (siz_wvl_pngau>0) then
1947        list_int(ii:ii+siz_wvl_pngau-1)=pawtab%wvl%pngau(1:siz_wvl_pngau)
1948        ii=ii+siz_wvl_pngau
1949      end if
1950      list_int(ii)=siz_wvl_rholoc_rad  ;ii=ii+1
1951      list_int(ii)=siz_wvl_rholoc_d  ;ii=ii+1
1952      list_int(ii)=pawtab%wvl%rholoc%msz  ;ii=ii+1
1953    end if
1954 
1955 !Then the data initialized later
1956 !...................................
1957    if (full_broadcast) then
1958 
1959 !Sizes of arrays
1960      list_int(ii)=siz_indklmn  ;ii=ii+1
1961      list_int(ii)=siz_klmntomn  ;ii=ii+1
1962      list_int(ii)=siz_kmix  ;ii=ii+1
1963      list_int(ii)=siz_lnproju  ;ii=ii+1
1964      list_int(ii)=siz_dltij  ;ii=ii+1
1965      list_int(ii)=siz_dshpfunc  ;ii=ii+1
1966      list_int(ii)=siz_eijkl  ;ii=ii+1
1967      list_int(ii)=siz_eijkl_sr  ;ii=ii+1
1968      list_int(ii)=siz_euijkl  ;ii=ii+1
1969      list_int(ii)=siz_euij_fll  ;ii=ii+1
1970      list_int(ii)=siz_fk  ;ii=ii+1
1971      list_int(ii)=siz_gammaij ;ii=ii+1
1972      list_int(ii)=siz_gnorm  ;ii=ii+1
1973      list_int(ii)=siz_nabla_ij  ;ii=ii+1
1974      list_int(ii)=siz_nabla_im_ij  ;ii=ii+1
1975      list_int(ii)=siz_nablaphi; ii=ii+1
1976      list_int(ii)=siz_phiphj  ;ii=ii+1
1977      list_int(ii)=siz_phiphjint  ;ii=ii+1
1978      list_int(ii)=siz_ph0phiint  ;ii=ii+1
1979      list_int(ii)=siz_qgrid_shp  ;ii=ii+1
1980      list_int(ii)=siz_qijl  ;ii=ii+1
1981      list_int(ii)=siz_rad_for_spline  ;ii=ii+1
1982      list_int(ii)=siz_shapefncg  ;ii=ii+1
1983      list_int(ii)=siz_sij  ;ii=ii+1
1984      list_int(ii)=siz_tnablaphi; ii=ii+1
1985      list_int(ii)=siz_tphitphj  ;ii=ii+1
1986      list_int(ii)=siz_vee  ;ii=ii+1
1987      list_int(ii)=siz_vex  ;ii=ii+1
1988      list_int(ii)=siz_zioneff  ;ii=ii+1
1989 !Integers
1990      list_int(ii)=pawtab%ij_proj  ;ii=ii+1
1991      list_int(ii)=pawtab%lcut_size  ;ii=ii+1
1992      list_int(ii)=pawtab%lexexch  ;ii=ii+1
1993      list_int(ii)=pawtab%lmnmix_sz  ;ii=ii+1
1994      list_int(ii)=pawtab%lpawu  ;ii=ii+1
1995      list_int(ii)=pawtab%mqgrid_shp  ;ii=ii+1
1996      list_int(ii)=pawtab%nproju  ;ii=ii+1
1997      list_int(ii)=pawtab%option_interaction_pawu ;ii=ii+1
1998      list_int(ii)=pawtab%useexexch  ;ii=ii+1
1999      list_int(ii)=pawtab%usepawu  ;ii=ii+1
2000      list_int(ii)=pawtab%usepotzero ;ii=ii+1
2001      list_int(ii)=pawtab%usespnorb ;ii=ii+1
2002 !Integer arrays
2003      if (siz_indklmn>0) then
2004        list_int(ii:ii+siz_indklmn-1)=reshape(pawtab%indklmn,(/siz_indklmn/))
2005        ii=ii+siz_indklmn
2006      end if
2007      if (siz_klmntomn>0) then
2008        list_int(ii:ii+siz_klmntomn-1)=reshape(pawtab%klmntomn,(/siz_klmntomn/))
2009        ii=ii+siz_klmntomn
2010      end if
2011      if (siz_kmix>0) then
2012        list_int(ii:ii+siz_kmix-1)=pawtab%kmix(1:siz_kmix)
2013        ii=ii+siz_kmix
2014      end if
2015      if (siz_lnproju>0) then
2016        list_int(ii:ii+siz_lnproju-1)=pawtab%lnproju(1:siz_lnproju)
2017        ii=ii+siz_lnproju
2018      end if
2019    end if ! full_broadcast
2020    ii=ii-1
2021 
2022    if (ii/=nn_int+nn_int_arr) then
2023      msg='the number of loaded integers is not correct!'
2024      LIBPAW_BUG(msg)
2025    end if
2026 
2027  end if ! me=0
2028 
2029 !Perfom the communication
2030 !-------------------------------------------------------------------------
2031 
2032  call xmpi_bcast(list_int,0,comm_mpi,ierr)
2033 
2034 !Fill the receiver from the buffer
2035 !-------------------------------------------------------------------------
2036  if (me/=0) then
2037    ii=1
2038 
2039 !First the data read from a psp file
2040 !...................................
2041 
2042 !Sizes of arrays (read from psp file)
2043    siz_indlmn=list_int(ii)  ;ii=ii+1
2044    siz_orbitals=list_int(ii)  ;ii=ii+1
2045    siz_coredens=list_int(ii)  ;ii=ii+1
2046    siz_coretau=list_int(ii)  ;ii=ii+1
2047    siz_dij0=list_int(ii)  ;ii=ii+1
2048    siz_kij=list_int(ii)  ;ii=ii+1
2049    siz_fock=list_int(ii)  ;ii=ii+1
2050    siz_phi=list_int(ii)  ;ii=ii+1
2051    siz_rhoij0=list_int(ii)  ;ii=ii+1
2052    siz_shape_alpha=list_int(ii)  ;ii=ii+1
2053    siz_shape_q=list_int(ii)  ;ii=ii+1
2054    siz_shapefunc=list_int(ii)  ;ii=ii+1
2055    siz_tcoredens=list_int(ii)  ;ii=ii+1
2056    siz_tcoretau=list_int(ii)  ;ii=ii+1
2057    siz_tcorespl=list_int(ii)  ;ii=ii+1
2058    siz_tcoretauspl=list_int(ii)  ;ii=ii+1
2059    siz_tphi=list_int(ii)  ;ii=ii+1
2060    siz_tproj=list_int(ii)  ;ii=ii+1
2061    siz_tvalespl=list_int(ii)  ;ii=ii+1
2062    siz_vhtnzc=list_int(ii)  ;ii=ii+1
2063    siz_vhnzc=list_int(ii)  ;ii=ii+1
2064    siz_vminushalf=list_int(ii)  ;ii=ii+1
2065    siz_wvlpaw=list_int(ii)  ;ii=ii+1
2066 !Integers (read from psp file)
2067    pawtab%basis_size=list_int(ii)  ;ii=ii+1
2068    pawtab%has_fock=list_int(ii)  ;ii=ii+1
2069    pawtab%has_kij=list_int(ii)  ;ii=ii+1
2070    pawtab%has_shapefncg=list_int(ii)  ;ii=ii+1
2071    pawtab%has_nabla=list_int(ii)  ;ii=ii+1
2072    pawtab%has_nablaphi=list_int(ii) ; ii=ii+1
2073    pawtab%has_tproj=list_int(ii)  ;ii=ii+1
2074    pawtab%has_tvale=list_int(ii)  ;ii=ii+1
2075    pawtab%has_coretau=list_int(ii)  ;ii=ii+1
2076    pawtab%has_vhtnzc=list_int(ii)  ;ii=ii+1
2077    pawtab%has_vhnzc=list_int(ii)  ;ii=ii+1
2078    pawtab%has_vminushalf=list_int(ii)  ;ii=ii+1
2079    pawtab%has_wvl=list_int(ii)  ;ii=ii+1
2080    pawtab%ij_size=list_int(ii)  ;ii=ii+1
2081    pawtab%l_size=list_int(ii)  ;ii=ii+1
2082    pawtab%lmn_size=list_int(ii)  ;ii=ii+1
2083    pawtab%lmn2_size=list_int(ii)  ;ii=ii+1
2084    pawtab%mesh_size=list_int(ii)  ;ii=ii+1
2085    pawtab%partialwave_mesh_size=list_int(ii)  ;ii=ii+1
2086    pawtab%core_mesh_size=list_int(ii)  ;ii=ii+1
2087    pawtab%coretau_mesh_size=list_int(ii)  ;ii=ii+1
2088    pawtab%vminus_mesh_size=list_int(ii)  ;ii=ii+1
2089    pawtab%tnvale_mesh_size=list_int(ii)  ;ii=ii+1
2090    pawtab%mqgrid=list_int(ii)  ;ii=ii+1
2091    pawtab%shape_lambda=list_int(ii)  ;ii=ii+1
2092    pawtab%shape_type=list_int(ii)  ;ii=ii+1
2093    pawtab%usetcore=list_int(ii)  ;ii=ii+1
2094    pawtab%usexcnhat=list_int(ii)  ;ii=ii+1
2095 !Integer arrays (read from psp file)
2096    if (allocated(pawtab%indlmn)) then
2097      LIBPAW_DEALLOCATE(pawtab%indlmn)
2098    end if
2099    if (siz_indlmn>0) then
2100      LIBPAW_ALLOCATE(pawtab%indlmn,(6,pawtab%lmn_size))
2101      pawtab%indlmn=reshape(list_int(ii:ii+siz_indlmn-1),(/6,pawtab%lmn_size/))
2102      ii=ii+siz_indlmn
2103    end if
2104    if (allocated(pawtab%orbitals)) then
2105      LIBPAW_DEALLOCATE(pawtab%orbitals)
2106    end if
2107    if (siz_orbitals>0) then
2108      LIBPAW_ALLOCATE(pawtab%orbitals,(pawtab%basis_size))
2109      pawtab%orbitals=list_int(ii:ii+pawtab%basis_size-1)
2110      ii=ii+siz_orbitals
2111    end if
2112 !Integers in datastructures (read from psp file)
2113    if (siz_wvlpaw==1) then
2114      call wvlpaw_allocate(pawtab%wvl)
2115      siz_wvl_pngau=list_int(ii)  ;ii=ii+1
2116      siz_wvl_parg=list_int(ii)  ;ii=ii+1
2117      siz_wvl_pfac=list_int(ii)  ;ii=ii+1
2118      pawtab%wvl%npspcode_init_guess=list_int(ii)  ;ii=ii+1
2119      pawtab%wvl%ptotgau=list_int(ii)  ;ii=ii+1
2120      if (allocated(pawtab%wvl%pngau)) then
2121        LIBPAW_DEALLOCATE(pawtab%wvl%pngau)
2122      end if
2123      if (siz_wvl_pngau>0) then
2124        LIBPAW_ALLOCATE(pawtab%wvl%pngau,(pawtab%basis_size))
2125        pawtab%wvl%pngau=list_int(ii:ii+pawtab%basis_size-1)
2126        ii=ii+siz_wvl_pngau
2127      end if
2128      siz_wvl_rholoc_rad=list_int(ii)  ;ii=ii+1
2129      siz_wvl_rholoc_d=list_int(ii)  ;ii=ii+1
2130      pawtab%wvl%rholoc%msz=list_int(ii)  ;ii=ii+1
2131    end if
2132 
2133 !Then the data initialized later
2134 !...................................
2135    if (full_broadcast) then
2136 
2137 !Sizes of arrays
2138      siz_indklmn=list_int(ii)  ;ii=ii+1
2139      siz_klmntomn=list_int(ii)  ;ii=ii+1
2140      siz_kmix=list_int(ii)  ;ii=ii+1
2141      siz_lnproju=list_int(ii)  ;ii=ii+1
2142      siz_dltij=list_int(ii)  ;ii=ii+1
2143      siz_dshpfunc=list_int(ii)  ;ii=ii+1
2144      siz_eijkl=list_int(ii)  ;ii=ii+1
2145      siz_eijkl_sr=list_int(ii)  ;ii=ii+1
2146      siz_euijkl=list_int(ii)  ;ii=ii+1
2147      siz_euij_fll=list_int(ii)  ;ii=ii+1
2148      siz_fk=list_int(ii)  ;ii=ii+1
2149      siz_gammaij=list_int(ii)  ;ii=ii+1
2150      siz_gnorm=list_int(ii)  ;ii=ii+1
2151      siz_nabla_ij=list_int(ii)  ;ii=ii+1
2152      siz_nabla_im_ij=list_int(ii)  ;ii=ii+1
2153      siz_nablaphi=list_int(ii)  ;ii=ii+1
2154      siz_phiphj=list_int(ii)  ;ii=ii+1
2155      siz_phiphjint=list_int(ii)  ;ii=ii+1
2156      siz_ph0phiint=list_int(ii)  ;ii=ii+1
2157      siz_qgrid_shp=list_int(ii)  ;ii=ii+1
2158      siz_qijl=list_int(ii)  ;ii=ii+1
2159      siz_rad_for_spline=list_int(ii)  ;ii=ii+1
2160      siz_shapefncg=list_int(ii)  ;ii=ii+1
2161      siz_sij=list_int(ii)  ;ii=ii+1
2162      siz_tnablaphi=list_int(ii)  ;ii=ii+1
2163      siz_tphitphj=list_int(ii)  ;ii=ii+1
2164      siz_vee=list_int(ii)  ;ii=ii+1
2165      siz_vex=list_int(ii)  ;ii=ii+1
2166      siz_zioneff=list_int(ii)  ;ii=ii+1
2167 !Integers
2168      pawtab%ij_proj=list_int(ii)  ;ii=ii+1
2169      pawtab%lcut_size=list_int(ii)  ;ii=ii+1
2170      pawtab%lexexch=list_int(ii)  ;ii=ii+1
2171      pawtab%lmnmix_sz=list_int(ii)  ;ii=ii+1
2172      pawtab%lpawu=list_int(ii)  ;ii=ii+1
2173      pawtab%mqgrid_shp=list_int(ii)  ;ii=ii+1
2174      pawtab%nproju=list_int(ii)  ;ii=ii+1
2175      pawtab%option_interaction_pawu=list_int(ii) ;ii=ii+1
2176      pawtab%useexexch=list_int(ii)  ;ii=ii+1
2177      pawtab%usepawu=list_int(ii)  ;ii=ii+1
2178      pawtab%usepotzero=list_int(ii) ;ii=ii+1
2179      pawtab%usespnorb=list_int(ii) ;ii=ii+1
2180 !Integer arrays
2181      if (allocated(pawtab%indklmn)) then
2182        LIBPAW_DEALLOCATE(pawtab%indklmn)
2183      end if
2184      if (siz_indklmn>0) then
2185        LIBPAW_ALLOCATE(pawtab%indklmn,(8,pawtab%lmn2_size))
2186        pawtab%indklmn=reshape(list_int(ii:ii+siz_indklmn-1),(/8,pawtab%lmn2_size/))
2187        ii=ii+siz_indklmn
2188      end if
2189      if (allocated(pawtab%klmntomn)) then
2190        LIBPAW_DEALLOCATE(pawtab%klmntomn)
2191      end if
2192      if (siz_klmntomn>0) then
2193        LIBPAW_ALLOCATE(pawtab%klmntomn,(4,pawtab%lmn2_size))
2194        pawtab%klmntomn=reshape(list_int(ii:ii+siz_klmntomn-1),(/4,pawtab%lmn2_size/))
2195        ii=ii+siz_klmntomn
2196      end if
2197      if (allocated(pawtab%kmix)) then
2198        LIBPAW_DEALLOCATE(pawtab%kmix)
2199      end if
2200      if (siz_kmix>0) then
2201        LIBPAW_ALLOCATE(pawtab%kmix,(pawtab%lmnmix_sz))
2202        pawtab%kmix=list_int(ii:ii+pawtab%lmnmix_sz-1)
2203        ii=ii+siz_kmix
2204      end if
2205      if (allocated(pawtab%lnproju)) then
2206        LIBPAW_DEALLOCATE(pawtab%lnproju)
2207      end if
2208      if (siz_lnproju>0) then
2209        LIBPAW_ALLOCATE(pawtab%lnproju,(pawtab%nproju))
2210        pawtab%lnproju=list_int(ii:ii+pawtab%nproju-1)
2211        ii=ii+siz_lnproju
2212      end if
2213    end if ! full_broadcast
2214    ii=ii-1
2215 
2216    if (ii/=nn_int+nn_int_arr) then
2217      msg='the number of broadcasted integers is not correct!'
2218      LIBPAW_BUG(msg)
2219    end if
2220 
2221  end if ! me/=0
2222  LIBPAW_DEALLOCATE(list_int)
2223 
2224 !Broadcast all the reals
2225 !=========================================================================
2226 
2227  LIBPAW_ALLOCATE(list_dpr,(nn_dpr+nn_dpr_arr))
2228 
2229 !Fill the buffer of the sender
2230 !-------------------------------------------------------------------------
2231  if (me==0) then
2232    ii=1
2233 
2234 !First the data read from a psp file
2235 !...................................
2236 
2237 !Reals (read from psp file)
2238    list_dpr(ii)=pawtab%beta    ;ii=ii+1
2239    list_dpr(ii)=pawtab%dncdq0  ;ii=ii+1
2240    list_dpr(ii)=pawtab%d2ncdq0  ;ii=ii+1
2241    list_dpr(ii)=pawtab%dnvdq0  ;ii=ii+1
2242    list_dpr(ii)=pawtab%dtaucdq0  ;ii=ii+1
2243    list_dpr(ii)=pawtab%ex_cc   ;ii=ii+1
2244    list_dpr(ii)=pawtab%exccore  ;ii=ii+1
2245    list_dpr(ii)=pawtab%lamb_shielding  ;ii=ii+1
2246    list_dpr(ii)=pawtab%rpaw  ;ii=ii+1
2247    list_dpr(ii)=pawtab%rshp  ;ii=ii+1
2248    list_dpr(ii)=pawtab%rcore  ;ii=ii+1
2249    list_dpr(ii)=pawtab%rcoretau  ;ii=ii+1
2250    list_dpr(ii)=pawtab%shape_sigma  ;ii=ii+1
2251 !Reals arrays (read from psp file)
2252    if (siz_coredens>0) then
2253      list_dpr(ii:ii+siz_coredens-1)=pawtab%coredens(1:siz_coredens)
2254      ii=ii+siz_coredens
2255    end if
2256    if (siz_coretau>0) then
2257      list_dpr(ii:ii+siz_coretau-1)=pawtab%coretau(1:siz_coretau)
2258      ii=ii+siz_coretau
2259    end if
2260    if (siz_dij0>0) then
2261      list_dpr(ii:ii+siz_dij0-1)=pawtab%dij0(1:siz_dij0)
2262      ii=ii+siz_dij0
2263    end if
2264    if (siz_fock>0) then
2265      list_dpr(ii:ii+siz_fock-1)=pawtab%ex_cvij(1:siz_fock)
2266      ii=ii+siz_fock
2267    end if
2268    if (siz_kij>0) then
2269      list_dpr(ii:ii+siz_kij-1)=pawtab%kij(1:siz_kij)
2270      ii=ii+siz_kij
2271    end if
2272    if (siz_phi>0) then
2273      list_dpr(ii:ii+siz_phi-1)=reshape(pawtab%phi,(/siz_phi/))
2274      ii=ii+siz_phi
2275    end if
2276    if (siz_rhoij0>0) then
2277      list_dpr(ii:ii+siz_rhoij0-1)=pawtab%rhoij0(1:siz_rhoij0)
2278      ii=ii+siz_rhoij0
2279    end if
2280    if (siz_shape_alpha>0) then
2281      list_dpr(ii:ii+siz_shape_alpha-1)=reshape(pawtab%shape_alpha,(/siz_shape_alpha/))
2282      ii=ii+siz_shape_alpha
2283    end if
2284    if (siz_shape_q>0) then
2285      list_dpr(ii:ii+siz_shape_q-1)=reshape(pawtab%shape_q,(/siz_shape_q/))
2286      ii=ii+siz_shape_q
2287    end if
2288    if (siz_shapefunc>0) then
2289      list_dpr(ii:ii+siz_shapefunc-1)=reshape(pawtab%shapefunc,(/siz_shapefunc/))
2290      ii=ii+siz_shapefunc
2291    end if
2292    if (siz_tcoredens>0) then
2293      list_dpr(ii:ii+siz_tcoredens-1)=reshape(pawtab%tcoredens,(/siz_tcoredens/))
2294      ii=ii+siz_tcoredens
2295    end if
2296    if (siz_tcoretau>0) then
2297      list_dpr(ii:ii+siz_tcoretau-1)=reshape(pawtab%tcoretau,(/siz_tcoretau/))
2298      ii=ii+siz_tcoretau
2299    end if
2300    if (siz_tcorespl>0) then
2301      list_dpr(ii:ii+siz_tcorespl-1)=reshape(pawtab%tcorespl,(/siz_tcorespl/))
2302      ii=ii+siz_tcorespl
2303    end if
2304    if (siz_tcoretauspl>0) then
2305      list_dpr(ii:ii+siz_tcoretauspl-1)=reshape(pawtab%tcoretauspl,(/siz_tcoretauspl/))
2306      ii=ii+siz_tcoretauspl
2307    end if
2308    if (siz_tphi>0) then
2309      list_dpr(ii:ii+siz_tphi-1)=reshape(pawtab%tphi,(/siz_tphi/))
2310      ii=ii+siz_tphi
2311    end if
2312    if (siz_tproj>0) then
2313      list_dpr(ii:ii+siz_tproj-1)=reshape(pawtab%tproj,(/siz_tproj/))
2314      ii=ii+siz_tproj
2315    end if
2316    if (siz_tvalespl>0) then
2317      list_dpr(ii:ii+siz_tvalespl-1)=reshape(pawtab%tvalespl,(/siz_tvalespl/))
2318      ii=ii+siz_tvalespl
2319    end if
2320    if (siz_vhtnzc>0) then
2321      list_dpr(ii:ii+siz_vhtnzc-1)=pawtab%vhtnzc(1:siz_vhtnzc)
2322      ii=ii+siz_vhtnzc
2323    end if
2324    if (siz_vhnzc>0) then
2325      list_dpr(ii:ii+siz_vhnzc-1)=pawtab%vhnzc(1:siz_vhnzc)
2326      ii=ii+siz_vhnzc
2327    end if
2328    if (siz_vminushalf>0) then
2329      list_dpr(ii:ii+siz_vminushalf-1)=pawtab%vminushalf(1:siz_vminushalf)
2330      ii=ii+siz_vminushalf
2331    end if
2332 !Reals in datastructures (read from psp file)
2333    if (siz_wvlpaw==1) then
2334      if (siz_wvl_parg>0) then
2335        list_dpr(ii:ii+siz_wvl_parg-1)=reshape(pawtab%wvl%parg,(/siz_wvl_parg/))
2336        ii=ii+siz_wvl_parg
2337      end if
2338      if (siz_wvl_pfac>0) then
2339        list_dpr(ii:ii+siz_wvl_pfac-1)=reshape(pawtab%wvl%pfac,(/siz_wvl_pfac/))
2340        ii=ii+siz_wvl_pfac
2341      end if
2342      if (siz_wvl_rholoc_rad>0) then
2343         list_dpr(ii:ii+siz_wvl_rholoc_rad-1)=pawtab%wvl%rholoc%rad(1:siz_wvl_rholoc_rad)
2344         ii=ii+siz_wvl_rholoc_rad
2345      end if
2346      if (siz_wvl_rholoc_d>0) then
2347         list_dpr(ii:ii+siz_wvl_rholoc_d-1)=reshape(pawtab%wvl%rholoc%d,(/siz_wvl_rholoc_d/))
2348         ii=ii+siz_wvl_rholoc_d
2349      end if
2350    end if
2351 
2352 !Then the data initialized later
2353 !...................................
2354    if (full_broadcast) then
2355 
2356 !Reals
2357      list_dpr(ii)=pawtab%exchmix  ;ii=ii+1
2358      list_dpr(ii)=pawtab%f4of2_sla  ;ii=ii+1
2359      list_dpr(ii)=pawtab%f6of2_sla  ;ii=ii+1
2360      list_dpr(ii)=pawtab%jpawu  ;ii=ii+1
2361      list_dpr(ii)=pawtab%upawu  ;ii=ii+1
2362 !Reals arrays
2363      if (siz_dltij>0) then
2364        list_dpr(ii:ii+siz_dltij-1)=pawtab%dltij(1:siz_dltij)
2365        ii=ii+siz_dltij
2366      end if
2367      if (siz_dshpfunc>0) then
2368        list_dpr(ii:ii+siz_dshpfunc-1)=reshape(pawtab%dshpfunc,(/siz_dshpfunc/))
2369        ii=ii+siz_dshpfunc
2370      end if
2371      if (siz_eijkl>0) then
2372        list_dpr(ii:ii+siz_eijkl-1)=reshape(pawtab%eijkl,(/siz_eijkl/))
2373        ii=ii+siz_eijkl
2374      end if
2375      if (siz_eijkl_sr>0) then
2376        list_dpr(ii:ii+siz_eijkl_sr-1)=reshape(pawtab%eijkl_sr,(/siz_eijkl_sr/))
2377        ii=ii+siz_eijkl_sr
2378      end if
2379      if (siz_euijkl>0) then
2380        list_dpr(ii:ii+siz_euijkl-1)=reshape(pawtab%euijkl,(/siz_euijkl/))
2381        ii=ii+siz_euijkl
2382      end if
2383      if (siz_euij_fll>0) then
2384        list_dpr(ii:ii+siz_euij_fll-1)=reshape(pawtab%euij_fll,(/siz_euij_fll/))
2385        ii=ii+siz_euij_fll
2386      end if
2387      if (siz_fk>0) then
2388        list_dpr(ii:ii+siz_fk-1)=reshape(pawtab%fk,(/siz_fk/))
2389        ii=ii+siz_fk
2390      end if
2391      if (siz_gammaij>0) then
2392        list_dpr(ii:ii+siz_gammaij-1)=pawtab%gammaij(1:siz_gammaij)
2393        ii=ii+siz_gammaij
2394      end if
2395      if (siz_gnorm>0) then
2396        list_dpr(ii:ii+siz_gnorm-1)=pawtab%gnorm(1:siz_gnorm)
2397        ii=ii+siz_gnorm
2398      end if
2399      if (siz_nabla_ij>0) then
2400        list_dpr(ii:ii+siz_nabla_ij-1)=reshape(pawtab%nabla_ij,(/siz_nabla_ij/))
2401        ii=ii+siz_nabla_ij
2402      end if
2403      if (siz_nabla_im_ij>0) then
2404        list_dpr(ii:ii+siz_nabla_im_ij-1)=reshape(pawtab%nabla_im_ij,(/siz_nabla_im_ij/))
2405        ii=ii+siz_nabla_im_ij
2406      end if
2407      if (siz_nablaphi>0) then
2408        list_dpr(ii:ii+siz_nablaphi-1)=reshape(pawtab%nablaphi,(/siz_nablaphi/))
2409        ii=ii+siz_nablaphi
2410      end if
2411      if (siz_phiphj>0) then
2412        list_dpr(ii:ii+siz_phiphj-1)=reshape(pawtab%phiphj,(/siz_phiphj/))
2413        ii=ii+siz_phiphj
2414      end if
2415      if (siz_phiphjint>0) then
2416        list_dpr(ii:ii+siz_phiphjint-1)=pawtab%phiphjint(1:siz_phiphjint)
2417        ii=ii+siz_phiphjint
2418      end if
2419      if (siz_ph0phiint>0) then
2420        list_dpr(ii:ii+siz_ph0phiint-1)=pawtab%ph0phiint(1:siz_ph0phiint)
2421        ii=ii+siz_ph0phiint
2422      end if
2423      if (siz_qgrid_shp>0) then
2424        list_dpr(ii:ii+siz_qgrid_shp-1)=pawtab%qgrid_shp(1:siz_qgrid_shp)
2425        ii=ii+siz_qgrid_shp
2426      end if
2427      if (siz_qijl>0) then
2428        list_dpr(ii:ii+siz_qijl-1)=reshape(pawtab%qijl,(/siz_qijl/))
2429        ii=ii+siz_qijl
2430      end if
2431      if (siz_rad_for_spline>0) then
2432        list_dpr(ii:ii+siz_rad_for_spline-1)=pawtab%rad_for_spline(1:siz_rad_for_spline)
2433        ii=ii+siz_rad_for_spline
2434      end if
2435      if (siz_shapefncg>0) then
2436        list_dpr(ii:ii+siz_shapefncg-1)=reshape(pawtab%shapefncg,(/siz_shapefncg/))
2437        ii=ii+siz_shapefncg
2438      end if
2439      if (siz_sij>0) then
2440        list_dpr(ii:ii+siz_sij-1)=pawtab%sij(1:siz_sij)
2441        ii=ii+siz_sij
2442      end if
2443      if (siz_tnablaphi>0) then
2444        list_dpr(ii:ii+siz_tnablaphi-1)=reshape(pawtab%tnablaphi,(/siz_tnablaphi/))
2445        ii=ii+siz_tnablaphi
2446      end if
2447      if (siz_tphitphj>0) then
2448        list_dpr(ii:ii+siz_tphitphj-1)=reshape(pawtab%tphitphj,(/siz_tphitphj/))
2449        ii=ii+siz_tphitphj
2450      end if
2451      if (siz_vee>0) then
2452        list_dpr(ii:ii+siz_vee-1)=reshape(pawtab%vee,(/siz_vee/))
2453        ii=ii+siz_vee
2454      end if
2455      if (siz_vex>0) then
2456        list_dpr(ii:ii+siz_vex-1)=reshape(pawtab%vex,(/siz_vex/))
2457        ii=ii+siz_vex
2458      end if
2459      if (siz_zioneff>0) then
2460        list_dpr(ii:ii+siz_zioneff-1)=pawtab%zioneff(1:siz_zioneff)
2461        ii=ii+siz_zioneff
2462      end if
2463 
2464    end if ! full_broadcast
2465    ii=ii-1
2466    if (ii/=nn_dpr+nn_dpr_arr) then
2467      msg='the number of loaded reals is not correct!'
2468      LIBPAW_BUG(msg)
2469    end if
2470 
2471  end if ! me=0
2472 
2473 !Perfom the communication
2474 !-------------------------------------------------------------------------
2475 
2476  call xmpi_bcast(list_dpr,0,comm_mpi,ierr)
2477 
2478 !Fill the receiver from the buffer
2479 !-------------------------------------------------------------------------
2480  if (me/=0) then
2481    ii=1
2482 
2483 !First the data read from a psp file
2484 !...................................
2485 
2486 !Reals (read from psp file)
2487    pawtab%beta=list_dpr(ii)    ;ii=ii+1
2488    pawtab%dncdq0=list_dpr(ii)  ;ii=ii+1
2489    pawtab%d2ncdq0=list_dpr(ii)  ;ii=ii+1
2490    pawtab%dnvdq0=list_dpr(ii)  ;ii=ii+1
2491    pawtab%dtaucdq0=list_dpr(ii)  ;ii=ii+1
2492    pawtab%ex_cc=list_dpr(ii)  ;ii=ii+1
2493    pawtab%exccore=list_dpr(ii)  ;ii=ii+1
2494    pawtab%lamb_shielding=list_dpr(ii)  ;ii=ii+1
2495    pawtab%rpaw=list_dpr(ii)  ;ii=ii+1
2496    pawtab%rshp=list_dpr(ii)  ;ii=ii+1
2497    pawtab%rcore=list_dpr(ii)  ;ii=ii+1
2498    pawtab%rcoretau=list_dpr(ii)  ;ii=ii+1
2499    pawtab%shape_sigma=list_dpr(ii)  ;ii=ii+1
2500 !Reals arrays (read from psp file)
2501    if (allocated(pawtab%coredens)) then
2502      LIBPAW_DEALLOCATE(pawtab%coredens)
2503    end if
2504    if (siz_coredens>0) then
2505      LIBPAW_ALLOCATE(pawtab%coredens,(pawtab%core_mesh_size))
2506      pawtab%coredens=list_dpr(ii:ii+pawtab%core_mesh_size-1)
2507      ii=ii+siz_coredens
2508    end if
2509    if (allocated(pawtab%coretau)) then
2510      LIBPAW_DEALLOCATE(pawtab%coretau)
2511    end if
2512    if (siz_coretau>0) then
2513      LIBPAW_ALLOCATE(pawtab%coretau,(pawtab%coretau_mesh_size))
2514      pawtab%coretau=list_dpr(ii:ii+pawtab%coretau_mesh_size-1)
2515      ii=ii+siz_coretau
2516    end if
2517    if (allocated(pawtab%dij0)) then
2518      LIBPAW_DEALLOCATE(pawtab%dij0)
2519    end if
2520    if (siz_dij0>0) then
2521      LIBPAW_ALLOCATE(pawtab%dij0,(pawtab%lmn2_size))
2522      pawtab%dij0=list_dpr(ii:ii+pawtab%lmn2_size-1)
2523      ii=ii+siz_dij0
2524    end if
2525    if (allocated(pawtab%ex_cvij)) then
2526      LIBPAW_DEALLOCATE(pawtab%ex_cvij)
2527    end if
2528    if (siz_fock>0) then
2529      LIBPAW_ALLOCATE(pawtab%ex_cvij,(pawtab%lmn2_size))
2530      pawtab%ex_cvij=list_dpr(ii:ii+pawtab%lmn2_size-1)
2531      ii=ii+siz_fock
2532    end if
2533    if (allocated(pawtab%kij)) then
2534      LIBPAW_DEALLOCATE(pawtab%kij)
2535    end if
2536    if (siz_kij>0) then
2537      LIBPAW_ALLOCATE(pawtab%kij,(pawtab%lmn2_size))
2538      pawtab%kij=list_dpr(ii:ii+pawtab%lmn2_size-1)
2539      ii=ii+siz_kij
2540    end if
2541    if (allocated(pawtab%phi)) then
2542      LIBPAW_DEALLOCATE(pawtab%phi)
2543    end if
2544    if (siz_phi>0) then
2545      LIBPAW_ALLOCATE(pawtab%phi,(pawtab%partialwave_mesh_size,pawtab%basis_size))
2546      pawtab%phi=reshape(list_dpr(ii:ii+siz_phi-1),(/pawtab%partialwave_mesh_size,pawtab%basis_size/))
2547      ii=ii+siz_phi
2548    end if
2549    if (allocated(pawtab%rhoij0)) then
2550      LIBPAW_DEALLOCATE(pawtab%rhoij0)
2551    end if
2552    if (siz_rhoij0>0) then
2553      LIBPAW_ALLOCATE(pawtab%rhoij0,(pawtab%lmn2_size))
2554      pawtab%rhoij0=list_dpr(ii:ii+pawtab%lmn2_size-1)
2555      ii=ii+siz_rhoij0
2556    end if
2557    if (allocated(pawtab%shape_alpha)) then
2558      LIBPAW_DEALLOCATE(pawtab%shape_alpha)
2559    end if
2560    if (siz_shape_alpha>0) then
2561      LIBPAW_ALLOCATE(pawtab%shape_alpha,(2,pawtab%l_size))
2562      pawtab%shape_alpha=reshape(list_dpr(ii:ii+siz_shape_alpha-1),(/2,pawtab%l_size/))
2563      ii=ii+siz_shape_alpha
2564    end if
2565    if (allocated(pawtab%shape_q)) then
2566      LIBPAW_DEALLOCATE(pawtab%shape_q)
2567    end if
2568    if (siz_shape_q>0) then
2569      LIBPAW_ALLOCATE(pawtab%shape_q,(2,pawtab%l_size))
2570      pawtab%shape_q=reshape(list_dpr(ii:ii+siz_shape_q-1),(/2,pawtab%l_size/))
2571      ii=ii+siz_shape_q
2572    end if
2573    if (allocated(pawtab%shapefunc)) then
2574      LIBPAW_DEALLOCATE(pawtab%shapefunc)
2575    end if
2576    if (siz_shapefunc>0) then
2577      LIBPAW_ALLOCATE(pawtab%shapefunc,(pawtab%mesh_size,pawtab%l_size))
2578      pawtab%shapefunc=reshape(list_dpr(ii:ii+siz_shapefunc-1),(/pawtab%mesh_size,pawtab%l_size/))
2579      ii=ii+siz_shapefunc
2580    end if
2581    if (allocated(pawtab%tcoredens)) then
2582      LIBPAW_DEALLOCATE(pawtab%tcoredens)
2583    end if
2584    if (siz_tcoredens>0) then
2585      sz2=siz_tcoredens/pawtab%core_mesh_size
2586      LIBPAW_ALLOCATE(pawtab%tcoredens,(pawtab%core_mesh_size,sz2))
2587      pawtab%tcoredens=reshape(list_dpr(ii:ii+siz_tcoredens-1),(/pawtab%core_mesh_size,sz2/))
2588      ii=ii+siz_tcoredens
2589    end if
2590    if (allocated(pawtab%tcoretau)) then
2591      LIBPAW_DEALLOCATE(pawtab%tcoretau)
2592    end if
2593    if (siz_tcoretau>0) then
2594      LIBPAW_ALLOCATE(pawtab%tcoretau,(pawtab%coretau_mesh_size))
2595      pawtab%tcoretau=list_dpr(ii:ii+siz_tcoretau-1)
2596      ii=ii+siz_tcoretau
2597    end if
2598    if (allocated(pawtab%tcorespl)) then
2599      LIBPAW_DEALLOCATE(pawtab%tcorespl)
2600    end if
2601    if (siz_tcorespl>0) then
2602      LIBPAW_ALLOCATE(pawtab%tcorespl,(pawtab%mqgrid,2))
2603      pawtab%tcorespl=reshape(list_dpr(ii:ii+siz_tcorespl-1),(/pawtab%mqgrid,2/))
2604      ii=ii+siz_tcorespl
2605    end if
2606    if (allocated(pawtab%tcoretauspl)) then
2607      LIBPAW_DEALLOCATE(pawtab%tcoretauspl)
2608    end if
2609    if (siz_tcoretauspl>0) then
2610      LIBPAW_ALLOCATE(pawtab%tcoretauspl,(pawtab%mqgrid,2))
2611      pawtab%tcoretauspl=reshape(list_dpr(ii:ii+siz_tcoretauspl-1),(/pawtab%mqgrid,2/))
2612      ii=ii+siz_tcoretauspl
2613    else
2614      LIBPAW_ALLOCATE(pawtab%tcoretauspl,(pawtab%mqgrid,0))
2615    end if
2616    if (allocated(pawtab%tphi)) then
2617      LIBPAW_DEALLOCATE(pawtab%tphi)
2618    end if
2619    if (siz_tphi>0) then
2620      LIBPAW_ALLOCATE(pawtab%tphi,(pawtab%partialwave_mesh_size,pawtab%basis_size))
2621      pawtab%tphi=reshape(list_dpr(ii:ii+siz_tphi-1),(/pawtab%partialwave_mesh_size,pawtab%basis_size/))
2622      ii=ii+siz_tphi
2623    end if
2624    if (allocated(pawtab%tproj)) then
2625      LIBPAW_DEALLOCATE(pawtab%tproj)
2626    end if
2627    if (siz_tproj>0) then
2628      sz1=siz_tproj/pawtab%basis_size
2629      LIBPAW_ALLOCATE(pawtab%tproj,(sz1,pawtab%basis_size))
2630      pawtab%tproj=reshape(list_dpr(ii:ii+siz_tproj-1),(/sz1,pawtab%basis_size/))
2631      ii=ii+siz_tproj
2632    end if
2633    if (allocated(pawtab%tvalespl)) then
2634      LIBPAW_DEALLOCATE(pawtab%tvalespl)
2635    end if
2636    if (siz_tvalespl>0) then
2637      sz1=siz_tvalespl/2
2638      LIBPAW_ALLOCATE(pawtab%tvalespl,(sz1,2))
2639      pawtab%tvalespl=reshape(list_dpr(ii:ii+siz_tvalespl-1),(/sz1,2/))
2640      ii=ii+siz_tvalespl
2641    end if
2642    if (allocated(pawtab%vhtnzc)) then
2643      LIBPAW_DEALLOCATE(pawtab%vhtnzc)
2644    end if
2645    if (siz_vhtnzc>0) then
2646      LIBPAW_ALLOCATE(pawtab%vhtnzc,(pawtab%mesh_size))
2647      pawtab%vhtnzc=list_dpr(ii:ii+pawtab%mesh_size-1)
2648      ii=ii+siz_vhtnzc
2649    end if
2650    if (allocated(pawtab%vhnzc)) then
2651      LIBPAW_DEALLOCATE(pawtab%vhnzc)
2652    end if
2653    if (siz_vhnzc>0) then
2654      LIBPAW_ALLOCATE(pawtab%vhnzc,(pawtab%mesh_size))
2655      pawtab%vhnzc=list_dpr(ii:ii+pawtab%mesh_size-1)
2656      ii=ii+siz_vhnzc
2657    end if
2658    if (allocated(pawtab%vminushalf)) then
2659      LIBPAW_DEALLOCATE(pawtab%vminushalf)
2660    end if
2661    if (siz_vminushalf>0) then
2662      LIBPAW_ALLOCATE(pawtab%vminushalf,(pawtab%mesh_size))
2663      pawtab%vminushalf=list_dpr(ii:ii+pawtab%mesh_size-1)
2664      ii=ii+siz_vminushalf
2665    end if
2666 !Reals in datastructures (read from psp file)
2667    if (siz_wvlpaw==1) then
2668      if (allocated(pawtab%wvl%parg)) then
2669        LIBPAW_DEALLOCATE(pawtab%wvl%parg)
2670      end if
2671      if (siz_wvl_parg>0) then
2672        LIBPAW_ALLOCATE(pawtab%wvl%parg,(2,pawtab%wvl%ptotgau))
2673        pawtab%wvl%parg=reshape(list_dpr(ii:ii+siz_wvl_parg-1),(/2,pawtab%wvl%ptotgau/))
2674        ii=ii+siz_wvl_parg
2675      end if
2676      if (allocated(pawtab%wvl%pfac)) then
2677        LIBPAW_DEALLOCATE(pawtab%wvl%pfac)
2678      end if
2679      if (siz_wvl_pfac>0) then
2680        LIBPAW_ALLOCATE(pawtab%wvl%pfac,(2,pawtab%wvl%ptotgau))
2681        pawtab%wvl%pfac=reshape(list_dpr(ii:ii+siz_wvl_pfac-1),(/2,pawtab%wvl%ptotgau/))
2682        ii=ii+siz_wvl_pfac
2683      end if
2684      if (allocated(pawtab%wvl%rholoc%rad)) then
2685         LIBPAW_DEALLOCATE(pawtab%wvl%rholoc%rad)
2686      end if
2687      if (siz_wvl_rholoc_rad>0) then
2688         sz1=pawtab%wvl%rholoc%msz
2689         LIBPAW_ALLOCATE(pawtab%wvl%rholoc%rad,(sz1))
2690         pawtab%wvl%rholoc%rad=list_dpr(ii:ii+sz1-1)
2691         ii=ii+siz_wvl_rholoc_rad
2692      end if
2693      if (allocated(pawtab%wvl%rholoc%d)) then
2694         LIBPAW_DEALLOCATE(pawtab%wvl%rholoc%d)
2695      end if
2696      if (siz_wvl_rholoc_d>0) then
2697         sz1=pawtab%wvl%rholoc%msz
2698         LIBPAW_ALLOCATE(pawtab%wvl%rholoc%d,(sz1,4))
2699         pawtab%wvl%rholoc%d=reshape(list_dpr(ii:ii+siz_wvl_rholoc_d-1),(/sz1,4/))
2700         ii=ii+siz_wvl_rholoc_d
2701      end if
2702    end if
2703 
2704 !Then the data initialized later
2705 !...................................
2706    if (full_broadcast) then
2707 
2708 !Reals
2709      pawtab%exchmix=list_dpr(ii)  ;ii=ii+1
2710      pawtab%f4of2_sla=list_dpr(ii)  ;ii=ii+1
2711      pawtab%f6of2_sla=list_dpr(ii)  ;ii=ii+1
2712      pawtab%jpawu=list_dpr(ii)  ;ii=ii+1
2713      pawtab%upawu=list_dpr(ii)  ;ii=ii+1
2714 !Reals arrays
2715      if (allocated(pawtab%dltij)) then
2716        LIBPAW_DEALLOCATE(pawtab%dltij)
2717      end if
2718      if (siz_dltij>0) then
2719        LIBPAW_ALLOCATE(pawtab%dltij,(pawtab%lmn2_size))
2720        pawtab%dltij=list_dpr(ii:ii+pawtab%lmn2_size-1)
2721        ii=ii+siz_dltij
2722      end if
2723      if (allocated(pawtab%dshpfunc)) then
2724        LIBPAW_DEALLOCATE(pawtab%dshpfunc)
2725      end if
2726      if (siz_dshpfunc>0) then
2727        LIBPAW_ALLOCATE(pawtab%dshpfunc,(pawtab%mesh_size,pawtab%l_size,4))
2728        pawtab%dshpfunc=reshape(list_dpr(ii:ii+siz_dshpfunc-1),(/pawtab%mesh_size,pawtab%l_size,4/))
2729        ii=ii+siz_dshpfunc
2730      end if
2731      if (allocated(pawtab%eijkl)) then
2732        LIBPAW_DEALLOCATE(pawtab%eijkl)
2733      end if
2734      if (siz_eijkl>0) then
2735        LIBPAW_ALLOCATE(pawtab%eijkl,(pawtab%lmn2_size,pawtab%lmn2_size))
2736        pawtab%eijkl=reshape(list_dpr(ii:ii+siz_eijkl-1),(/pawtab%lmn2_size,pawtab%lmn2_size/))
2737        ii=ii+siz_eijkl
2738      end if
2739      if (allocated(pawtab%eijkl_sr)) then
2740        LIBPAW_DEALLOCATE(pawtab%eijkl_sr)
2741      end if
2742      if (siz_eijkl_sr>0) then
2743        LIBPAW_ALLOCATE(pawtab%eijkl_sr,(pawtab%lmn2_size,pawtab%lmn2_size))
2744        pawtab%eijkl_sr=reshape(list_dpr(ii:ii+siz_eijkl_sr-1),(/pawtab%lmn2_size,pawtab%lmn2_size/))
2745        ii=ii+siz_eijkl_sr
2746      end if
2747      if (allocated(pawtab%euijkl)) then
2748        LIBPAW_DEALLOCATE(pawtab%euijkl)
2749      end if
2750      if (siz_euijkl>0) then
2751        LIBPAW_ALLOCATE(pawtab%euijkl,(3,pawtab%lmn_size,pawtab%lmn_size,pawtab%lmn_size,pawtab%lmn_size))
2752        pawtab%euijkl=reshape(list_dpr(ii:ii+siz_euijkl-1),(/3,pawtab%lmn_size,pawtab%lmn_size,pawtab%lmn_size,pawtab%lmn_size/))
2753        ii=ii+siz_euijkl
2754      end if
2755      if (allocated(pawtab%euij_fll)) then
2756        LIBPAW_DEALLOCATE(pawtab%euij_fll)
2757      end if
2758      if (siz_euij_fll>0) then
2759        LIBPAW_ALLOCATE(pawtab%euij_fll,(pawtab%lmn2_size))
2760        pawtab%euij_fll=reshape(list_dpr(ii:ii+siz_euij_fll-1),(/pawtab%lmn2_size/))
2761        ii=ii+siz_euij_fll
2762      end if
2763      if (allocated(pawtab%fk)) then
2764        LIBPAW_DEALLOCATE(pawtab%fk)
2765      end if
2766      if (siz_fk>0) then
2767        LIBPAW_ALLOCATE(pawtab%fk,(6,4))
2768        pawtab%fk=reshape(list_dpr(ii:ii+siz_fk-1),(/6,4/))
2769        ii=ii+siz_fk
2770      end if
2771      if (allocated(pawtab%gammaij)) then
2772        LIBPAW_DEALLOCATE(pawtab%gammaij)
2773      end if
2774      if (siz_gammaij>0) then
2775        LIBPAW_ALLOCATE(pawtab%gammaij,(pawtab%l_size))
2776        pawtab%gammaij=list_dpr(ii:ii+pawtab%l_size-1)
2777        ii=ii+siz_gammaij
2778      end if
2779      if (allocated(pawtab%gnorm)) then
2780        LIBPAW_DEALLOCATE(pawtab%gnorm)
2781      end if
2782      if (siz_gnorm>0) then
2783        LIBPAW_ALLOCATE(pawtab%gnorm,(pawtab%l_size))
2784        pawtab%gnorm=list_dpr(ii:ii+pawtab%l_size-1)
2785        ii=ii+siz_gnorm
2786      end if
2787      if (allocated(pawtab%nabla_ij)) then
2788        LIBPAW_DEALLOCATE(pawtab%nabla_ij)
2789      end if
2790      if (siz_nabla_ij>0) then
2791        LIBPAW_ALLOCATE(pawtab%nabla_ij,(3,pawtab%lmn_size,pawtab%lmn_size))
2792        pawtab%nabla_ij=reshape(list_dpr(ii:ii+siz_nabla_ij-1),(/3,pawtab%lmn_size,pawtab%lmn_size/))
2793        ii=ii+siz_nabla_ij
2794      end if
2795      if (allocated(pawtab%nabla_im_ij)) then
2796        LIBPAW_DEALLOCATE(pawtab%nabla_im_ij)
2797      end if
2798      if (siz_nabla_im_ij>0) then
2799        LIBPAW_ALLOCATE(pawtab%nabla_im_ij,(3,pawtab%lmn_size,pawtab%lmn_size))
2800        pawtab%nabla_im_ij=reshape(list_dpr(ii:ii+siz_nabla_im_ij-1),(/3,pawtab%lmn_size,pawtab%lmn_size/))
2801        ii=ii+siz_nabla_im_ij
2802      end if
2803      if (allocated(pawtab%nablaphi)) then
2804        LIBPAW_DEALLOCATE(pawtab%nablaphi)
2805      end if
2806      if (siz_nablaphi>0) then
2807        LIBPAW_ALLOCATE(pawtab%nablaphi,(pawtab%partialwave_mesh_size,pawtab%basis_size))
2808        pawtab%nablaphi=reshape(list_dpr(ii:ii+siz_nablaphi-1),(/pawtab%partialwave_mesh_size,pawtab%basis_size/))
2809        ii=ii+siz_nablaphi
2810      end if
2811      if (allocated(pawtab%phiphj)) then
2812        LIBPAW_DEALLOCATE(pawtab%phiphj)
2813      end if
2814      if (siz_phiphj>0) then
2815        LIBPAW_ALLOCATE(pawtab%phiphj,(pawtab%mesh_size,pawtab%ij_size))
2816        pawtab%phiphj=reshape(list_dpr(ii:ii+siz_phiphj-1),(/pawtab%mesh_size,pawtab%ij_size/))
2817        ii=ii+siz_phiphj
2818      end if
2819      if (allocated(pawtab%phiphjint)) then
2820        LIBPAW_DEALLOCATE(pawtab%phiphjint)
2821      end if
2822      if (siz_phiphjint>0) then
2823        LIBPAW_ALLOCATE(pawtab%phiphjint,(pawtab%ij_proj))
2824        pawtab%phiphjint=list_dpr(ii:ii+pawtab%ij_proj-1)
2825        ii=ii+siz_phiphjint
2826      end if
2827      if (allocated(pawtab%ph0phiint)) then
2828        LIBPAW_DEALLOCATE(pawtab%ph0phiint)
2829      end if
2830      if (siz_ph0phiint>0) then
2831        LIBPAW_ALLOCATE(pawtab%ph0phiint,(pawtab%ij_proj))
2832        pawtab%ph0phiint=list_dpr(ii:ii+pawtab%ij_proj-1)
2833        ii=ii+siz_ph0phiint
2834      end if
2835      if (allocated(pawtab%qgrid_shp)) then
2836        LIBPAW_DEALLOCATE(pawtab%qgrid_shp)
2837      end if
2838      if (siz_qgrid_shp>0) then
2839        LIBPAW_ALLOCATE(pawtab%qgrid_shp,(pawtab%mqgrid_shp))
2840        pawtab%qgrid_shp=list_dpr(ii:ii+pawtab%mqgrid_shp-1)
2841        ii=ii+siz_qgrid_shp
2842      end if
2843      if (allocated(pawtab%qijl)) then
2844        LIBPAW_DEALLOCATE(pawtab%qijl)
2845      end if
2846      if (siz_qijl>0) then
2847        LIBPAW_ALLOCATE(pawtab%qijl,(pawtab%l_size**2,pawtab%lmn2_size))
2848        pawtab%qijl=reshape(list_dpr(ii:ii+siz_qijl-1),(/pawtab%l_size**2,pawtab%lmn2_size/))
2849        ii=ii+siz_qijl
2850      end if
2851      if (allocated(pawtab%rad_for_spline)) then
2852        LIBPAW_DEALLOCATE(pawtab%rad_for_spline)
2853      end if
2854      if (siz_rad_for_spline>0) then
2855        LIBPAW_ALLOCATE(pawtab%rad_for_spline,(pawtab%mesh_size))
2856        pawtab%rad_for_spline=list_dpr(ii:ii+pawtab%mesh_size-1)
2857        ii=ii+siz_rad_for_spline
2858      end if
2859      if (allocated(pawtab%shapefncg)) then
2860        LIBPAW_DEALLOCATE(pawtab%shapefncg)
2861      end if
2862      if (siz_shapefncg>0) then
2863        LIBPAW_ALLOCATE(pawtab%shapefncg,(pawtab%mqgrid_shp,2,pawtab%l_size))
2864        pawtab%shapefncg=reshape(list_dpr(ii:ii+siz_shapefncg-1),(/pawtab%mqgrid_shp,2,pawtab%l_size/))
2865        ii=ii+siz_shapefncg
2866      end if
2867      if (allocated(pawtab%sij)) then
2868        LIBPAW_DEALLOCATE(pawtab%sij)
2869      end if
2870      if (siz_sij>0) then
2871        LIBPAW_ALLOCATE(pawtab%sij,(pawtab%lmn2_size))
2872        pawtab%sij=list_dpr(ii:ii+pawtab%lmn2_size-1)
2873        ii=ii+siz_sij
2874      end if
2875     if (allocated(pawtab%tnablaphi)) then
2876        LIBPAW_DEALLOCATE(pawtab%tnablaphi)
2877      end if
2878      if (siz_tnablaphi>0) then
2879        LIBPAW_ALLOCATE(pawtab%tnablaphi,(pawtab%partialwave_mesh_size,pawtab%basis_size))
2880        pawtab%tphi=reshape(list_dpr(ii:ii+siz_tnablaphi-1),(/pawtab%partialwave_mesh_size,pawtab%basis_size/))
2881        ii=ii+siz_tnablaphi
2882      end if
2883      if (allocated(pawtab%tphitphj)) then
2884        LIBPAW_DEALLOCATE(pawtab%tphitphj)
2885      end if
2886      if (siz_tphitphj>0) then
2887        LIBPAW_ALLOCATE(pawtab%tphitphj,(pawtab%mesh_size,pawtab%ij_size))
2888        pawtab%tphitphj=reshape(list_dpr(ii:ii+siz_tphitphj-1),(/pawtab%mesh_size,pawtab%ij_size/))
2889        ii=ii+siz_tphitphj
2890      end if
2891      if (allocated(pawtab%vee)) then
2892        LIBPAW_DEALLOCATE(pawtab%vee)
2893      end if
2894      if (siz_vee>0) then
2895        sz1=2*pawtab%lpawu+1
2896        LIBPAW_ALLOCATE(pawtab%vee,(sz1,sz1,sz1,sz1))
2897        pawtab%vee=reshape(list_dpr(ii:ii+siz_vee-1),(/sz1,sz1,sz1,sz1/))
2898        ii=ii+siz_vee
2899      end if
2900      if (allocated(pawtab%vex)) then
2901        LIBPAW_DEALLOCATE(pawtab%vex)
2902      end if
2903      if (siz_vex>0) then
2904        sz1=2*pawtab%lexexch+1
2905        LIBPAW_ALLOCATE(pawtab%vex,(sz1,sz1,sz1,sz1,4))
2906        pawtab%vex=reshape(list_dpr(ii:ii+siz_vex-1),(/sz1,sz1,sz1,sz1,4/))
2907        ii=ii+siz_vex
2908      end if
2909      if (allocated(pawtab%zioneff)) then
2910        LIBPAW_DEALLOCATE(pawtab%zioneff)
2911      end if
2912      if (siz_zioneff>0) then
2913        LIBPAW_ALLOCATE(pawtab%zioneff,(pawtab%ij_proj))
2914        pawtab%zioneff=list_dpr(ii:ii+pawtab%ij_proj-1)
2915        ii=ii+siz_zioneff
2916      end if
2917 
2918    end if ! full_broadcast
2919    ii=ii-1
2920 
2921    if (ii/=nn_dpr+nn_dpr_arr) then
2922      msg='the number of broadcasted reals is not correct!'
2923      LIBPAW_BUG(msg)
2924    end if
2925 
2926  end if ! me/=0
2927  LIBPAW_DEALLOCATE(list_dpr)
2928 
2929 end subroutine pawtab_bcast

m_pawtab/pawtab_free_0D [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  pawtab_free_0D

FUNCTION

  Deallocate pointers and nullify flags in a pawtab structure

SIDE EFFECTS

  Pawtab<type(pawtab_type)>=PAW arrays tabulated.
  All allocated arrays in Pawtab are deallocated

SOURCE

768 subroutine pawtab_free_0D(Pawtab)
769 
770 !Arguments ------------------------------------
771 !arrays
772  type(Pawtab_type),intent(inout) :: Pawtab
773 
774 !Local variables-------------------------------
775 
776 ! *************************************************************************
777 
778  !@Pawtab_type
779 
780  if (allocated(Pawtab%indklmn))  then
781   LIBPAW_DEALLOCATE(Pawtab%indklmn)
782  end if
783  if (allocated(Pawtab%indlmn))  then
784   LIBPAW_DEALLOCATE(Pawtab%indlmn)
785  end if
786  if (allocated(Pawtab%klmntomn))  then
787    LIBPAW_DEALLOCATE(Pawtab%klmntomn)
788  end if
789  if (allocated(Pawtab%kmix))  then
790    LIBPAW_DEALLOCATE(Pawtab%kmix)
791  end if
792  if (allocated(Pawtab%lnproju))  then
793    LIBPAW_DEALLOCATE(Pawtab%lnproju)
794  end if
795  if (allocated(Pawtab%coredens))  then
796    LIBPAW_DEALLOCATE(Pawtab%coredens)
797  end if
798  if (allocated(Pawtab%coretau))  then
799    LIBPAW_DEALLOCATE(Pawtab%coretau)
800  end if
801  if (allocated(Pawtab%dij0))  then
802    LIBPAW_DEALLOCATE(Pawtab%dij0)
803  end if
804  if (allocated(Pawtab%dltij))  then
805    LIBPAW_DEALLOCATE(Pawtab%dltij)
806  end if
807  if (allocated(Pawtab%dshpfunc))  then
808    LIBPAW_DEALLOCATE(Pawtab%dshpfunc)
809  end if
810  if (allocated(Pawtab%eijkl))  then
811    LIBPAW_DEALLOCATE(Pawtab%eijkl)
812  end if
813  if (allocated(Pawtab%eijkl_sr))  then
814    LIBPAW_DEALLOCATE(Pawtab%eijkl_sr)
815  end if
816  if (allocated(Pawtab%euijkl))  then
817    LIBPAW_DEALLOCATE(Pawtab%euijkl)
818  end if
819  if (allocated(Pawtab%euij_fll))  then
820    LIBPAW_DEALLOCATE(Pawtab%euij_fll)
821  end if
822  if (allocated(Pawtab%fk))  then
823    LIBPAW_DEALLOCATE(Pawtab%fk)
824  end if
825  if (allocated(Pawtab%gammaij))  then
826    LIBPAW_DEALLOCATE(Pawtab%gammaij)
827  end if
828  if (allocated(Pawtab%gnorm))  then
829    LIBPAW_DEALLOCATE(Pawtab%gnorm)
830  end if
831  if (allocated(Pawtab%ex_cvij))  then
832    LIBPAW_DEALLOCATE(Pawtab%ex_cvij)
833  end if
834  if (allocated(Pawtab%kij))  then
835    LIBPAW_DEALLOCATE(Pawtab%kij)
836  end if
837  if (allocated(Pawtab%nabla_ij))  then
838    LIBPAW_DEALLOCATE(Pawtab%nabla_ij)
839  end if
840  if (allocated(Pawtab%nabla_im_ij))  then
841    LIBPAW_DEALLOCATE(Pawtab%nabla_im_ij)
842  end if
843  if (allocated(Pawtab%nablaphi))  then
844    LIBPAW_DEALLOCATE(Pawtab%nablaphi)
845  end if
846  if (allocated(Pawtab%orbitals)) then
847    LIBPAW_DEALLOCATE(Pawtab%orbitals)
848  end if
849  if (allocated(Pawtab%phi))  then
850    LIBPAW_DEALLOCATE(Pawtab%phi)
851  end if
852  if (allocated(Pawtab%phiphj))  then
853    LIBPAW_DEALLOCATE(Pawtab%phiphj)
854  end if
855  if (allocated(Pawtab%phiphjint))  then
856    LIBPAW_DEALLOCATE(Pawtab%phiphjint)
857  end if
858  if (allocated(Pawtab%ph0phiint))  then
859    LIBPAW_DEALLOCATE(Pawtab%ph0phiint)
860  end if
861  if (allocated(Pawtab%qgrid_shp))  then
862    LIBPAW_DEALLOCATE(Pawtab%qgrid_shp)
863  end if
864  if (allocated(Pawtab%qijl))  then
865    LIBPAW_DEALLOCATE(Pawtab%qijl)
866  end if
867  if (allocated(Pawtab%rad_for_spline))  then
868    LIBPAW_DEALLOCATE(Pawtab%rad_for_spline)
869  end if
870  if (allocated(Pawtab%rhoij0))  then
871    LIBPAW_DEALLOCATE(Pawtab%rhoij0)
872  end if
873  if (allocated(Pawtab%shape_alpha))  then
874    LIBPAW_DEALLOCATE(Pawtab%shape_alpha)
875  end if
876  if (allocated(Pawtab%shape_q))  then
877    LIBPAW_DEALLOCATE(Pawtab%shape_q)
878  end if
879  if (allocated(Pawtab%shapefunc))  then
880    LIBPAW_DEALLOCATE(Pawtab%shapefunc)
881  end if
882  if (allocated(Pawtab%shapefncg))  then
883    LIBPAW_DEALLOCATE(Pawtab%shapefncg)
884  end if
885  if (allocated(Pawtab%sij))  then
886    LIBPAW_DEALLOCATE(Pawtab%sij)
887  end if
888  if (allocated(Pawtab%tcoredens))  then
889    LIBPAW_DEALLOCATE(Pawtab%tcoredens)
890  end if
891  if (allocated(Pawtab%tcoretau))  then
892    LIBPAW_DEALLOCATE(Pawtab%tcoretau)
893  end if
894  if (allocated(Pawtab%tcorespl))  then
895    LIBPAW_DEALLOCATE(Pawtab%tcorespl)
896  end if
897  if (allocated(Pawtab%tcoretauspl))  then
898    LIBPAW_DEALLOCATE(Pawtab%tcoretauspl)
899  end if
900  if (allocated(Pawtab%tnablaphi))  then
901    LIBPAW_DEALLOCATE(Pawtab%tnablaphi)
902  end if
903  if (allocated(Pawtab%tphi))  then
904    LIBPAW_DEALLOCATE(Pawtab%tphi)
905  end if
906  if (allocated(Pawtab%tphitphj))  then
907    LIBPAW_DEALLOCATE(Pawtab%tphitphj)
908  end if
909  if (allocated(Pawtab%tproj)) then
910    LIBPAW_DEALLOCATE(Pawtab%tproj)
911  end if
912  if (allocated(Pawtab%tvalespl))  then
913    LIBPAW_DEALLOCATE(Pawtab%tvalespl)
914  end if
915  if (allocated(Pawtab%vee))  then
916    LIBPAW_DEALLOCATE(Pawtab%vee)
917  end if
918  if (allocated(Pawtab%Vex))  then
919    LIBPAW_DEALLOCATE(Pawtab%Vex)
920  end if
921  if (allocated(Pawtab%vhtnzc))  then
922    LIBPAW_DEALLOCATE(Pawtab%vhtnzc)
923  end if
924  if (allocated(Pawtab%VHnZC))  then
925    LIBPAW_DEALLOCATE(Pawtab%VHnZC)
926  end if
927  if (allocated(Pawtab%vminushalf))  then
928    LIBPAW_DEALLOCATE(Pawtab%vminushalf)
929  end if
930  if (allocated(Pawtab%zioneff))  then
931    LIBPAW_DEALLOCATE(Pawtab%zioneff)
932  end if
933 
934  call wvlpaw_free(Pawtab%wvl)
935 
936  ! === Reset all flags and sizes ===
937 
938 !CAUTION: do not reset these flags
939 !They are set from input data and must be kept
940 !Pawtab%has_kij=0
941 !Pawtab%has_tproj=0
942 !Pawtab%has_coretau=0
943 !Pawtab%has_tvale=0
944 !Pawtab%has_vhtnzc=0
945 !Pawtab%has_vhnzc=0
946 !Pawtab%has_vminushalf=0
947 !Pawtab%has_nabla=0
948 !Pawtab%has_nablaphi=0
949 !Pawtab%has_shapefncg=0
950 !Pawtab%has_wvl=0
951 
952  Pawtab%usetcore=0
953  Pawtab%usexcnhat=0
954  Pawtab%useexexch=0
955  Pawtab%usepawu=0
956  Pawtab%usepotzero=0
957  Pawtab%usespnorb=0
958  Pawtab%mqgrid=0
959  Pawtab%mqgrid_shp=0
960 
961  Pawtab%basis_size=0
962  Pawtab%ij_proj=0
963  Pawtab%ij_size=0
964  Pawtab%lcut_size=0
965  Pawtab%l_size=0
966  Pawtab%lexexch=-1
967  Pawtab%lmn_size=0
968  Pawtab%lmn2_size=0
969  Pawtab%lmnmix_sz=0
970  Pawtab%lpawu=-1
971  Pawtab%nproju=0
972  Pawtab%option_interaction_pawu=0
973  Pawtab%mesh_size=0
974  Pawtab%partialwave_mesh_size=0
975  Pawtab%core_mesh_size=0
976  Pawtab%coretau_mesh_size=0
977  Pawtab%vminus_mesh_size=0
978  Pawtab%tnvale_mesh_size=0
979  Pawtab%shape_type=-10
980 
981 end subroutine pawtab_free_0D

m_pawtab/pawtab_free_1D [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  pawtab_free_1D

FUNCTION

  Destroy (deallocate) all pointers in an array of pawtab data structures

SOURCE

 995 subroutine pawtab_free_1D(Pawtab)
 996 
 997 !Arguments ------------------------------------
 998  type(pawtab_type),intent(inout) :: Pawtab(:)
 999 
1000 !Local variables-------------------------------
1001  integer :: ii,nn
1002 
1003 ! *************************************************************************
1004 
1005  !@pawtab_type
1006 
1007  nn=size(Pawtab)
1008  if (nn==0) return
1009 
1010  do ii=1,nn
1011    call pawtab_free_0D(Pawtab(ii))
1012  end do
1013 
1014 end subroutine pawtab_free_1D

m_pawtab/pawtab_nullify_0D [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  pawtab_nullify_0D

FUNCTION

  Nullify pointers and flags in a pawtab structure

SIDE EFFECTS

  Pawtab<type(pawtab_type)>=PAW arrays tabulated.
                            Nullified in output

SOURCE

659 subroutine pawtab_nullify_0D(Pawtab)
660 
661 !Arguments ------------------------------------
662 !arrays
663  type(Pawtab_type),intent(inout) :: Pawtab
664 
665 !Local variables-------------------------------
666 
667 ! *************************************************************************
668 
669  !@Pawtab_type
670  nullify(Pawtab%wvl)
671 
672  ! === Reset all flags and sizes ===
673 
674 !Flags controlling optional arrays
675  Pawtab%has_fock=0
676  Pawtab%has_kij=0
677  Pawtab%has_tproj=0
678  Pawtab%has_tvale=0
679  Pawtab%has_coretau=0
680  Pawtab%has_vhtnzc=0
681  Pawtab%has_vhnzc=0
682  Pawtab%has_vminushalf=0
683  Pawtab%has_nabla=0
684  Pawtab%has_nablaphi=0
685  Pawtab%has_shapefncg=0
686  Pawtab%has_wvl=0
687 
688  Pawtab%usetcore=0
689  Pawtab%usexcnhat=0
690  Pawtab%useexexch=0
691  Pawtab%usepawu=0
692  Pawtab%usepotzero=0
693  Pawtab%usespnorb=0
694  Pawtab%mqgrid=0
695  Pawtab%mqgrid_shp=0
696 
697  Pawtab%basis_size=0
698  Pawtab%ij_proj=0
699  Pawtab%ij_size=0
700  Pawtab%lcut_size=0
701  Pawtab%l_size=0
702  Pawtab%lexexch=-1
703  Pawtab%lmn_size=0
704  Pawtab%lmn2_size=0
705  Pawtab%lmnmix_sz=0
706  Pawtab%lpawu=-1
707  Pawtab%nproju=0
708  Pawtab%option_interaction_pawu=0
709  Pawtab%mesh_size=0
710  Pawtab%partialwave_mesh_size=0
711  Pawtab%core_mesh_size=0
712  Pawtab%coretau_mesh_size=0
713  Pawtab%vminus_mesh_size=0
714  Pawtab%tnvale_mesh_size=0
715  Pawtab%shape_type=-10
716 
717 end subroutine pawtab_nullify_0D

m_pawtab/pawtab_nullify_1D [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  pawtab_nullify_1D

FUNCTION

  Nullify all pointers in an array of pawtab data structures

SOURCE

731 subroutine pawtab_nullify_1D(Pawtab)
732 
733 !Arguments ------------------------------------
734  type(pawtab_type),intent(inout) :: Pawtab(:)
735 
736 !Local variables-------------------------------
737  integer :: ii,nn
738 
739 ! *************************************************************************
740 
741  !@pawtab_type
742 
743  nn=size(Pawtab)
744  if (nn==0) return
745 
746  do ii=1,nn
747    call pawtab_nullify_0D(Pawtab(ii))
748  end do
749 
750 end subroutine pawtab_nullify_1D

m_pawtab/pawtab_print [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

 pawtab_print

FUNCTION

  Print out the content of a pawtab datastructure

INPUTS

  Pawtab<pawtab_type> Only for PAW, TABulated data initialized at start

OUTPUT

  Only writing

SOURCE

1148 subroutine pawtab_print(Pawtab,header,unit,prtvol,mode_paral)
1149 
1150 !Arguments ------------------------------------
1151 !scalars
1152  integer,optional,intent(in) :: unit,prtvol
1153  character(len=4),optional,intent(in) :: mode_paral
1154  character(len=*),optional,intent(in) :: header
1155 !arrays
1156  type(Pawtab_type) :: Pawtab(:)
1157 
1158 !Local variables-------------------------------
1159 !scalars
1160  integer :: ityp,ntypat,my_unt,my_prtvol
1161  character(len=4) :: my_mode
1162  character(len=500) :: msg
1163 
1164 ! *************************************************************************
1165 
1166  my_unt   =ab_out ; if (PRESENT(unit      )) my_unt   =unit
1167  my_prtvol=0      ; if (PRESENT(prtvol    )) my_prtvol=prtvol
1168  my_mode  ='COLL' ; if (PRESENT(mode_paral)) my_mode  =mode_paral
1169 
1170  write(msg,'(6a)')&
1171 &  ' ==================================== ',ch10,&
1172 &  ' ==== Info on PAW TABulated data ==== ',ch10,&
1173 &  ' ==================================== ',ch10
1174  if (PRESENT(header)) msg=' ==== '//TRIM(ADJUSTL(header))//' ==== '
1175  call wrtout(my_unt,msg,my_mode)
1176 
1177  ntypat=SIZE(Pawtab(:))
1178 
1179  do ityp=1,ntypat
1180 
1181   ! Print out integer values (dimensions)
1182   write(msg,'(a)')'                                 '
1183   call wrtout(my_unt,msg,my_mode)
1184   write(msg,'(a)')'  ****************************** '
1185   call wrtout(my_unt,msg,my_mode)
1186   write(msg,'(a,i4,a)')'  **** Atom type ',ityp,' ****   '
1187   call wrtout(my_unt,msg,my_mode)
1188   write(msg,'(a)')'  ****************************** '
1189   call wrtout(my_unt,msg,my_mode)
1190   write(msg,'(a,i4)')'  Number of (n,l) elements ....................... ',Pawtab(ityp)%basis_size
1191   call wrtout(my_unt,msg,my_mode)
1192   write(msg,'(a,i4)')'  Number of (l,m,n) elements ..................... ',Pawtab(ityp)%lmn_size
1193   call wrtout(my_unt,msg,my_mode)
1194   write(msg,'(a,i4)')'  Number of (i,j) elements (packed form) ......... ',Pawtab(ityp)%ij_size
1195   call wrtout(my_unt,msg,my_mode)
1196   write(msg,'(a,i4)')'  Max L+1 leading to non-zero Gaunt .............. ',Pawtab(ityp)%l_size
1197   call wrtout(my_unt,msg,my_mode)
1198   write(msg,'(a,i4)')'  Max L+1 leading to non-zero Gaunt (pawlcutd) ... ',Pawtab(ityp)%lcut_size
1199   call wrtout(my_unt,msg,my_mode)
1200   write(msg,'(a,i4)')'  lmn2_size ...................................... ',Pawtab(ityp)%lmn2_size
1201   call wrtout(my_unt,msg,my_mode)
1202   write(msg,'(a,i4)')'  lmnmix_sz ...................................... ',Pawtab(ityp)%lmnmix_sz
1203   call wrtout(my_unt,msg,my_mode)
1204   write(msg,'(a,i4)')'  Size of radial mesh ............................ ',Pawtab(ityp)%mesh_size
1205   call wrtout(my_unt,msg,my_mode)
1206   write(msg,'(a,i4)')'  Size of radial mesh for partial waves........... ',Pawtab(ityp)%partialwave_mesh_size
1207   call wrtout(my_unt,msg,my_mode)
1208   write(msg,'(a,i4)')'  Size of radial mesh for [pseudo] core density... ',Pawtab(ityp)%core_mesh_size
1209   call wrtout(my_unt,msg,my_mode)
1210   write(msg,'(a,i4)')'  Size of radial mesh for [pseudo] kin core density',Pawtab(ityp)%coretau_mesh_size
1211   call wrtout(my_unt,msg,my_mode)
1212   write(msg,'(a,i4)')'  Size of radial mesh for pseudo valence density.. ',Pawtab(ityp)%tnvale_mesh_size
1213   call wrtout(my_unt,msg,my_mode)
1214   write(msg,'(a,i4)')'  No of Q-points for tcorespl/tvalespl/tcoretauspl ',Pawtab(ityp)%mqgrid
1215   call wrtout(my_unt,msg,my_mode)
1216   write(msg,'(a,i4)')'  No of Q-points for the radial shape functions .. ',Pawtab(ityp)%mqgrid_shp
1217   call wrtout(my_unt,msg,my_mode)
1218   write(msg,'(a,i4)')'  Radial shape function type ..................... ',Pawtab(ityp)%shape_type
1219   call wrtout(my_unt,msg,my_mode)
1220   write(msg,'(a,i4)')'  shape_lambda ................................... ',Pawtab(ityp)%shape_lambda
1221   call wrtout(my_unt,msg,my_mode)
1222   write(msg,'(a,i4)')'  Use pseudized core density ..................... ',Pawtab(ityp)%usetcore
1223   call wrtout(my_unt,msg,my_mode)
1224   write(msg,'(a,i4)')'  Option for the use of hat density in XC terms .. ',Pawtab(ityp)%usexcnhat
1225   call wrtout(my_unt,msg,my_mode)
1226   write(msg,'(a,i4)')'  Use DFT+U ...................................... ',Pawtab(ityp)%usepawu
1227   call wrtout(my_unt,msg,my_mode)
1228   if (Pawtab(ityp)%usepawu/=0) then
1229     write(msg,'(a,i4)')'  L on which U is applied ........................ ',Pawtab(ityp)%lpawu
1230     call wrtout(my_unt,msg,my_mode)
1231   end if
1232   write(msg,'(a,i4)')'  Use Local Exact exchange ....................... ',Pawtab(ityp)%useexexch
1233   call wrtout(my_unt,msg,my_mode)
1234   if (Pawtab(ityp)%useexexch/=0) then
1235     write(msg,'(a,i4)')'  L on which local exact-exchange is applied ..... ',Pawtab(ityp)%lexexch
1236     call wrtout(my_unt,msg,my_mode)
1237   end if
1238   if (Pawtab(ityp)%usepawu/=0.or.Pawtab(ityp)%useexexch/=0) then
1239     write(msg,'(a,i4)')'  Number of (i,j) elements for PAW+U or EXX ..... ',Pawtab(ityp)%ij_proj
1240     call wrtout(my_unt,msg,my_mode)
1241     write(msg,'(a,i4)')'  Number of projectors on which U or EXX acts .... ',Pawtab(ityp)%nproju
1242     call wrtout(my_unt,msg,my_mode)
1243     write(msg,'(a,i4)')'  Option interaction for PAW+U (double-counting).. ',Pawtab(ityp)%option_interaction_pawu
1244     call wrtout(my_unt,msg,my_mode)
1245   end if
1246   write(msg,'(a,i4)')'  Use potential zero ............................. ',Pawtab(ityp)%usepotzero
1247   call wrtout(my_unt,msg,my_mode)
1248   write(msg,'(a,i4)')'  Use spin-orbit coupling ........................ ',Pawtab(ityp)%usespnorb
1249   call wrtout(my_unt,msg,my_mode)
1250 
1251   ! "Has" flags
1252   write(msg,'(a,i4)')'  Has Fock  ...................................... ',Pawtab(ityp)%has_fock
1253   call wrtout(my_unt,msg,my_mode)
1254   write(msg,'(a,i4)')'  Has kij   ...................................... ',Pawtab(ityp)%has_kij
1255   call wrtout(my_unt,msg,my_mode)
1256   write(msg,'(a,i4)')'  Has tproj ...................................... ',Pawtab(ityp)%has_tproj
1257   call wrtout(my_unt,msg,my_mode)
1258   write(msg,'(a,i4)')'  Has tvale ...................................... ',Pawtab(ityp)%has_tvale
1259   call wrtout(my_unt,msg,my_mode)
1260   write(msg,'(a,i4)')'  Has coretau .................................... ',Pawtab(ityp)%has_coretau
1261   call wrtout(my_unt,msg,my_mode)
1262   write(msg,'(a,i4)')'  Has vhtnzc ..................................... ',Pawtab(ityp)%has_vhtnzc
1263   call wrtout(my_unt,msg,my_mode)
1264   write(msg,'(a,i4)')'  Has vhnzc ...................................... ',Pawtab(ityp)%has_vhnzc
1265   call wrtout(my_unt,msg,my_mode)
1266   write(msg,'(a,i4)')'  Has vminushalf ................................. ',Pawtab(ityp)%has_vminushalf
1267   call wrtout(my_unt,msg,my_mode)
1268   write(msg,'(a,i4)')'  Has nabla ...................................... ',Pawtab(ityp)%has_nabla
1269   call wrtout(my_unt,msg,my_mode)
1270   write(msg,'(a,i4)')'  Has nablaphi ................................... ',Pawtab(ityp)%has_nablaphi
1271   call wrtout(my_unt,msg,my_mode)
1272   write(msg,'(a,i4)')'  Has shapefuncg ................................. ',Pawtab(ityp)%has_shapefncg
1273   call wrtout(my_unt,msg,my_mode)
1274   write(msg,'(a,i4)')'  Has wvl ........................................ ',Pawtab(ityp)%has_wvl
1275   call wrtout(my_unt,msg,my_mode)
1276   !
1277   ! Real scalars
1278   write(msg,'(a,es16.8)')'  beta ............................................',Pawtab(ityp)%beta
1279   call wrtout(my_unt,msg,my_mode)
1280   write(msg,'(a,es16.8)')'  1/q d(tNcore(q))/dq for q=0 .....................',Pawtab(ityp)%dncdq0
1281   call wrtout(my_unt,msg,my_mode)
1282   write(msg,'(a,es16.8)')'  d^2(tNcore(q))/dq^2 for q=0 .....................',Pawtab(ityp)%d2ncdq0
1283   call wrtout(my_unt,msg,my_mode)
1284   write(msg,'(a,es16.8)')'  1/q d(tNvale(q))/dq for q=0 .....................',Pawtab(ityp)%dnvdq0
1285   call wrtout(my_unt,msg,my_mode)
1286   if (Pawtab(ityp)%has_coretau/=0) then
1287     write(msg,'(a,es16.8)')'  1/q d(tTAUcore(q))/dq for q=0 ...................',Pawtab(ityp)%dtaucdq0
1288     call wrtout(my_unt,msg,my_mode)
1289   end if
1290   if (Pawtab(ityp)%has_fock/=0) then
1291     write(msg,'(a,es16.8)')'  Core-core Fock energy  ..........................',Pawtab(ityp)%ex_cc
1292     call wrtout(my_unt,msg,my_mode)
1293   end if
1294   write(msg,'(a,es16.8)')'  XC energy for the core density ..................',Pawtab(ityp)%exccore
1295   call wrtout(my_unt,msg,my_mode)
1296   write(msg,'(a,es16.8)')'  Lamb shielding due to core density ..............',Pawtab(ityp)%lamb_shielding
1297   call wrtout(my_unt,msg,my_mode)
1298   write(msg,'(a,es16.8)')'  Radius of the PAW sphere ........................',Pawtab(ityp)%rpaw
1299   call wrtout(my_unt,msg,my_mode)
1300   write(msg,'(a,es16.8)')'  Compensation charge radius (if >rshp, g(r)=0) ...',Pawtab(ityp)%rshp !(if r>rshp, g(r)=zero)
1301   call wrtout(my_unt,msg,my_mode)
1302   if (Pawtab(ityp)%shape_type==2) then
1303    write(msg,'(a,es16.8)')'  Sigma parameter in gaussian shape function ......',Pawtab(ityp)%shape_sigma !(shape_type=2)
1304    call wrtout(my_unt,msg,my_mode)
1305   end if
1306   if (Pawtab(ityp)%usepawu/=0) then
1307    write(msg,'(a,es16.8)')'  Value of the U parameter [eV] ...................',Pawtab(ityp)%upawu*Ha_eV
1308    call wrtout(my_unt,msg,my_mode)
1309    write(msg,'(a,es16.8)')'  Value of the J parameter [eV] ...................',Pawtab(ityp)%jpawu*Ha_eV
1310    call wrtout(my_unt,msg,my_mode)
1311   end if
1312   if (Pawtab(ityp)%useexexch/=0) then
1313     write(msg,'(a,es16.8)')'  Mixing of exact exchange (PBE0) .................',Pawtab(ityp)%exchmix
1314     call wrtout(my_unt,msg,my_mode)
1315   end if
1316  if (associated(Pawtab(ityp)%wvl)) then
1317    write(msg,'(a,es16.8)')'  WARNING: This Pawtab structure contains WVL data.'
1318    call wrtout(my_unt,msg,my_mode)
1319  end if
1320 
1321  end do ! ityp
1322 
1323  ! The other (huge) arrays are not reported..
1324 
1325 end subroutine pawtab_print

m_pawtab/pawtab_set_flags_0D [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  pawtab_set_flags_0D

FUNCTION

  Set flags controlling optional arrays in a pawtab datastructure

SOURCE

1028 subroutine pawtab_set_flags_0D(Pawtab,has_coretau,has_fock,has_kij,has_tproj,has_tvale,has_vhnzc,&
1029 &                              has_vhtnzc,has_nabla,has_nablaphi,has_shapefncg,has_vminushalf,has_wvl)
1030 
1031 !Arguments ------------------------------------
1032  integer,intent(in),optional :: has_coretau,has_fock,has_kij,has_tproj,has_tvale
1033  integer,intent(in),optional :: has_vhnzc,has_vhtnzc,has_vminushalf
1034  integer,intent(in),optional :: has_nabla,has_nablaphi,has_shapefncg,has_wvl
1035  type(pawtab_type),intent(inout) :: Pawtab
1036 
1037 !Local variables-------------------------------
1038 
1039 ! *************************************************************************
1040 
1041  !@pawtab_type
1042 
1043  Pawtab%has_fock      =0
1044  Pawtab%has_kij       =0
1045  Pawtab%has_tproj     =0
1046  Pawtab%has_tvale     =0
1047  Pawtab%has_coretau   =0
1048  Pawtab%has_vhnzc     =0
1049  Pawtab%has_vhtnzc    =0
1050  Pawtab%has_nabla     =0
1051  Pawtab%has_nablaphi  =0
1052  Pawtab%has_shapefncg =0
1053  Pawtab%has_vminushalf=0
1054  Pawtab%has_wvl       =0
1055  if (present(has_fock))      Pawtab%has_fock=has_fock
1056  if (present(has_kij))       Pawtab%has_kij=has_kij
1057  if (present(has_tproj))     Pawtab%has_tproj=has_tproj
1058  if (present(has_tvale))     Pawtab%has_tvale=has_tvale
1059  if (present(has_coretau))   Pawtab%has_coretau=has_coretau
1060  if (present(has_vhnzc))     Pawtab%has_vhnzc=has_vhnzc
1061  if (present(has_vhtnzc))    Pawtab%has_vhtnzc=has_vhtnzc
1062  if (present(has_nabla))     Pawtab%has_nabla=has_nabla
1063  if (present(has_nablaphi))  Pawtab%has_nablaphi=has_nablaphi
1064  if (present(has_shapefncg) )Pawtab%has_shapefncg=has_shapefncg
1065  if (present(has_vminushalf))Pawtab%has_vminushalf=has_vminushalf
1066  if (present(has_wvl))       Pawtab%has_wvl=has_wvl
1067 
1068 end subroutine pawtab_set_flags_0D

m_pawtab/pawtab_set_flags_1D [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  pawtab_set_flags_1D

FUNCTION

  Set flags controlling optional arrays in an array of pawtab datastructures
 if (present(has_tvale))    Pawtab%has_tvale=has_tvale

SOURCE

1083 subroutine pawtab_set_flags_1D(Pawtab,has_coretau,has_fock,has_kij,has_tproj,has_tvale,has_vhnzc,&
1084 &                              has_vhtnzc,has_nabla,has_nablaphi,has_shapefncg,has_vminushalf,has_wvl)
1085 
1086 !Arguments ------------------------------------
1087  integer,intent(in),optional :: has_coretau,has_fock,has_kij,has_tproj,has_tvale,has_vhnzc,has_vhtnzc
1088  integer,intent(in),optional :: has_nabla,has_nablaphi,has_shapefncg,has_vminushalf,has_wvl
1089  type(pawtab_type),intent(inout) :: Pawtab(:)
1090 
1091 !Local variables-------------------------------
1092  integer :: ii,nn
1093 
1094 ! *************************************************************************
1095 
1096  !@pawtab_type
1097 
1098  nn=size(Pawtab)
1099  if (nn==0) return
1100 
1101  do ii=1,nn
1102    Pawtab(ii)%has_fock      =0
1103    Pawtab(ii)%has_kij       =0
1104    Pawtab(ii)%has_tproj     =0
1105    Pawtab(ii)%has_tvale     =0
1106    Pawtab(ii)%has_coretau   =0
1107    Pawtab(ii)%has_vhnzc     =0
1108    Pawtab(ii)%has_vhtnzc    =0
1109    Pawtab(ii)%has_nabla     =0
1110    Pawtab(ii)%has_nablaphi  =0
1111    Pawtab(ii)%has_shapefncg =0
1112    Pawtab(ii)%has_vminushalf=0
1113    Pawtab(ii)%has_wvl       =0
1114    if (present(has_fock))      Pawtab(ii)%has_fock=has_fock
1115    if (present(has_kij))       Pawtab(ii)%has_kij=has_kij
1116    if (present(has_tproj))     Pawtab(ii)%has_tproj=has_tproj
1117    if (present(has_tvale))     Pawtab(ii)%has_tvale=has_tvale
1118    if (present(has_coretau))   Pawtab(ii)%has_coretau=has_coretau
1119    if (present(has_vhnzc))     Pawtab(ii)%has_vhnzc=has_vhnzc
1120    if (present(has_vhtnzc))    Pawtab(ii)%has_vhtnzc=has_vhtnzc
1121    if (present(has_nabla))     Pawtab(ii)%has_nabla=has_nabla
1122    if (present(has_nablaphi))  Pawtab(ii)%has_nablaphi=has_nablaphi
1123    if (present(has_shapefncg)) Pawtab(ii)%has_shapefncg=has_shapefncg
1124    if (present(has_vminushalf))Pawtab(ii)%has_vminushalf=has_vminushalf
1125    if (present(has_wvl))       Pawtab(ii)%has_wvl=has_wvl
1126  end do
1127 
1128 end subroutine pawtab_set_flags_1D

m_pawtab/pawtab_type [ Types ]

[ Top ] [ m_pawtab ] [ Types ]

NAME

 pawtab_type

FUNCTION

 This structured datatype contains TABulated data for PAW (from pseudopotential)
 used in PAW calculations.

SOURCE

129  type,public :: pawtab_type
130 
131 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
132 ! declared in another part of ABINIT, that might need to take into account your modification.
133 
134 !Integer scalars
135 
136   integer :: basis_size
137    ! Number of elements for the paw nl basis on the considered atom type
138 
139   integer :: has_coretau
140    ! Flag controling use of core kinetic enrgy density (AE and pseudo)
141    ! if 1, [t]coretau() is allocated.
142    ! if 2, [t]coretau() is computed and stored.
143 
144   integer :: has_fock
145    ! if 1, onsite matrix elements of the core-valence Fock operator are allocated
146    ! if 2, onsite matrix elements of the core-valence Fock operator are computed and stored
147 
148   integer :: has_kij
149    ! if 1, onsite matrix elements of the kinetic operator are allocated
150    ! if 2, onsite matrix elements of the kinetic operator are computed and stored
151 
152   integer :: has_shapefncg
153    ! if 1, the spherical Fourier transforms of the radial shape functions are allocated
154    ! if 2, the spherical Fourier transforms of the radial shape functions are computed and stored
155 
156   integer :: has_nabla
157    ! if 1, onsite matrix elements of the nabla operator are allocated.
158    ! if 2, onsite matrix elements of the nabla operator are computed and stored.
159    ! if 3, onsite matrix elements of the nabla operator between valence and core states are computed and stored.
160    ! if 4, onsite matrix elements of the nabla operator between valence and core spinor states are are computed and stored.
161 
162   integer :: has_nablaphi
163    ! if 1, nablaphi are allocated
164    ! if 2, nablaphi are computed for MetaGGA and stored.
165 
166   integer :: has_tproj
167    ! Flag controling use of projectors in real space (0 if tnval is unknown)
168    ! if 1, tproj() is allocated.
169    ! if 2, tproj() is computed and stored.
170 
171   integer :: has_tvale
172    ! Flag controling use of pseudized valence density (0 if tnval is unknown)
173    ! if 1, tvalespl() is allocated.
174    ! if 2, tvalespl() is computed and stored.
175 
176   integer :: has_vhtnzc
177    ! if 1, space for vhtnzc is allocated
178    ! if 2, vhtnzc has been read from PAW file and stored
179 
180   integer :: has_vhnzc
181    ! if 1, space for vhnzc is allocated
182    ! if 2, vhnzc has been computed and stored
183 
184   integer :: has_vminushalf
185    ! has_vminushalf=0 ; vminushal is not allocated
186    ! has_vminushalf=1 ; vminushal is not allocated and stored
187 
188   integer :: has_wvl
189    ! if 1, data for wavelets (pawwvl) are allocated
190    ! if 2, data for wavelets (pawwvl) are computed and stored
191 
192   integer :: ij_proj
193    ! Number of (i,j) elements for the orbitals on which U acts (PAW+U only)
194    ! on the considered atom type (ij_proj=1 (1 projector), 3 (2 projectors)...)
195    ! Also used for local exact-exchange
196 
197   integer :: ij_size
198    ! Number of (i,j) elements for the symetric paw basis
199    ! on the considered atom type (ij_size=basis_size*(basis_size+1)/2)
200 
201   integer :: lcut_size
202    ! Maximum value of l+1 leading to non zero Gaunt coeffs
203    ! modified by dtset%pawlcutd
204    ! lcut_size=min(2*l_max,dtset%pawlcutd)+1
205 
206   integer :: l_size
207    ! Maximum value of l+1 leading to non zero Gaunt coeffs
208    ! l_size=2*l_max-1
209 
210   integer :: lexexch
211    ! lexexch gives l on which local exact-exchange is applied for a given type of atom.
212 
213   integer :: lmn_size
214    ! Number of (l,m,n) elements for the paw basis
215 
216   integer :: lmn2_size
217    ! lmn2_size=lmn_size*(lmn_size+1)/2
218    ! where lmn_size is the number of (l,m,n) elements for the paw basis
219 
220   integer :: lmnmix_sz
221    ! lmnmix_sz=number of klmn=(lmn,lmn_prime) verifying l<=lmix and l_prime<=lmix
222 
223   integer :: lpawu
224    ! lpawu gives l on which U is applied for a given type of atom.
225 
226   integer :: nproju
227    ! nproju is the number of projectors for orbitals on which paw+u acts.
228    ! Also used for local exact-exchange
229 
230   integer :: option_interaction_pawu
231    ! Option for interaction energy (PAW+U) in case of non-collinear magnetism:
232    ! 1: E_int=-J/4.N.(N-2)
233    ! 2: E_int=-J/2.(Nup.(Nup-1)+Ndn.(Ndn-1))    (Nup and Ndn are ill-defined)
234    ! 3: E_int=-J/4.( N.(N-2) + mx^2 + my^2 + mz^2 )
235 
236   integer :: mesh_size
237    ! Dimension of radial mesh for generic arrays contained in this pawtab datastructure
238    ! The mesh is usually defined up to the PAW augmentation region boundary
239    ! (+ a few additional points). May be different from pawrad%mesh_size
240 
241   integer :: core_mesh_size
242    ! Dimension of radial mesh for core density
243 
244   integer :: coretau_mesh_size
245    ! Dimension of radial mesh for core density
246 
247   integer :: vminus_mesh_size
248    ! Dimension of radial mesh for vminushalf
249 
250   integer :: partialwave_mesh_size
251    ! Dimension of radial mesh for partial waves (phi, tphi)
252    ! May be different from pawrad%mesh_size and pawtab%mesh_size
253 
254   integer :: tnvale_mesh_size
255    ! Dimension of radial mesh for tnvale
256 
257   integer :: mqgrid
258    ! Number of points in the reciprocal space grid on which
259    ! the radial functions (tcorespl, tvalespl...) are specified
260    ! Same as psps%mqgrid_vl
261 
262   integer :: mqgrid_shp
263    ! Number of points in the reciprocal space grid on which
264    ! the radial shape functions (shapefncg) are given
265 
266   integer :: usespnorb
267    ! usespnorb=0 ; no spin-orbit coupling
268    ! usespnorb=1 ; spin-orbit coupling
269 
270   integer :: shape_lambda
271    ! Lambda parameter in gaussian shapefunction (shape_type=2)
272 
273   integer :: shape_type
274    ! Radial shape function type
275    ! shape_type=-1 ; g(r)=numeric (read from psp file)
276    ! shape_type= 1 ; g(r)=[sin(pi*r/rshp)/(pi*r/rshp)]**2 if r<=rshp, zero if r>rshp
277    ! shape_type= 2 ; g(r)=exp[-(r/sigma)**lambda]
278    ! shape_type= 3 ; gl(r)=Alpha(1,l)*jl(q(1,l)*r)+Alpha(2,l)*jl(q(2,l)*r) for each l
279 
280   integer :: useexexch
281    ! useexexch=0 ; do not use local exact-exchange
282    ! useexexch=1 ; use local exact-exchange
283 
284   integer :: usepawu
285    ! usepawu= 0 ; do not use PAW+U formalism
286    ! usepawu= 1 ; use PAW+U formalism (Full localized limit)
287    ! usepawu= 2 ; use PAW+U formalism (Around Mean Field)
288    ! usepawu= 3 ; use PAW+U formalism (Around Mean Field) - Alternative
289    ! usepawu= 4 ; use PAW+U formalism (FLL) without polarization in the XC
290    ! usepawu=-1 ; use PAW+U formalism (FLL) - No use of the occupation matrix - Experimental
291    ! usepawu=-2 ; use PAW+U formalism (AMF) - No use of the occupation matrix - Experimental
292    ! usepawu=-4 ; use PAW+U formalism (FLL) without polarization in the XC - No use of the occupation matrix - Experimental
293    ! usepawu=10 ; use PAW+U within DMFT
294    ! usepawu=14 ; use PAW+U within DMFT without polarization in the XC
295 
296   integer :: usepotzero
297    ! usepotzero=0 if it is the Kresse-Joubert convention
298    ! usepotzero=1 if it is the new convention
299    ! usepotzero=2 if it is the PWscf convention
300 
301   integer :: usetcore
302    ! Flag controling use of pseudized core density (0 if tncore=zero)
303 
304   integer :: usexcnhat
305    ! 0 if compensation charge density is not included in XC terms
306    ! 1 if compensation charge density is included in XC terms
307 
308 !Real (real(dp)) scalars
309 
310   real(dp) :: beta
311    ! contains the integral of the difference between vH[nZc] and vH[tnZc]
312 
313   real(dp) :: dncdq0
314    ! Gives 1/q d(tNcore(q))/dq for q=0
315    ! (tNcore(q) = FT of pseudo core density)
316 
317   real(dp) :: d2ncdq0
318    ! Gives contribution of d2(tNcore(q))/d2q for q=0
319    ! \int{(16/15)*pi^5*n(r)*r^6* dr}
320    ! (tNcore(q) = FT of pseudo core density)
321 
322   real(dp) :: dnvdq0
323    ! Gives 1/q d(tNvale(q))/dq for q=0
324    ! (tNvale(q) = FT of pseudo valence density)
325 
326   real(dp) :: dtaucdq0
327    ! Gives 1/q d(tTAUcore(q))/dq for q=0
328    ! (tTAUcore(q) = FT of pseudo core kinetic density)
329 
330   real(dp) :: ex_cc
331    ! Exchange energy for the core-core interaction of the Fock operator
332 
333   real(dp) :: exccore
334    ! Exchange-correlation energy for the core density
335 
336   real(dp) :: exchmix
337    ! mixing of exact exchange; default is 0.25 (PBE0)
338 
339   real(dp) :: f4of2_sla
340    ! Ratio of Slater Integrals F4 and F2
341 
342   real(dp) :: f6of2_sla
343    ! Ratio of Slater Integrals F6 and F4
344 
345   real(dp) :: jpawu
346    ! Value of J parameter for paw+u for a given type.
347 
348    real(dp) :: lamb_shielding=0.0D0
349    ! Lamb shielding used in NMR shielding calcs (see m_orbmag.F90)
350 
351   real(dp) :: rpaw
352    ! Radius of PAW sphere
353 
354   real(dp) :: rshp
355    ! Compensation charge radius (if r>rshp, g(r)=zero)
356 
357   real(dp) :: rcore
358    ! Radius of core corrections (rcore >= rpaw)
359 
360   real(dp) :: rcoretau
361    ! Radius of kinetic core corrections (rcoretau >= rpaw)
362 
363   real(dp) :: shape_sigma
364    ! Sigma parameter in gaussian shapefunction (shape_type=2)
365 
366   real(dp) :: upawu
367    ! Value of U parameter for paw+u for a given type.
368 
369 !Objects
370   type(wvlpaw_type), pointer :: wvl
371    !variable containing objects needed
372    !for wvl+paw implementation
373    !Warning: it is a pointer; it has to be allocated before use
374 
375 !Integer arrays
376 
377   integer, allocatable :: indklmn(:,:)
378    ! indklmn(8,lmn2_size)
379    ! Array giving klm, kln, abs(il-jl), (il+jl), ilm and jlm, ilmn and jlmn for each klmn=(ilmn,jlmn)
380    ! Note: ilmn=(il,im,in) and ilmn<=jlmn
381 
382   integer, allocatable :: indlmn(:,:)
383    ! indlmn(6,lmn_size)
384    ! For each type of psp,
385    ! array giving l,m,n,lm,ln,spin for i=lmn (if useylm=1)
386 
387   integer, allocatable :: klmntomn(:,:)
388    ! klmntomn(4,lmn2_size)
389    ! Array giving im, jm ,in, and jn for each klmn=(ilmn,jlmn)
390    ! Note: ilmn=(il,im,in) and ilmn<=jlmn
391    ! NB: klmntomn is an application and not a bijection
392 
393   integer, allocatable :: kmix(:)
394    ! kmix(lmnmix_sz)
395    ! Indirect array selecting the klmn=(lmn,lmn_prime) verifying l<=lmix and l_prime<=lmix
396 
397   integer, allocatable :: lnproju(:)
398    ! lnproju(nproju) gives ln (index for phi) for each projectors on which U acts (PAW+U only)
399    ! nproju is 1 or 2 and  is the number of projectors for correlated orbitals
400    ! Also used for local exact-exchange
401 
402   integer, allocatable :: orbitals(:)
403    ! (basis_size)
404    ! gives the l quantum number per basis element
405 
406 !Real (real(dp)) arrays
407 
408   real(dp), allocatable :: coredens(:)
409    ! coredens(mesh_size)
410    ! Gives the core density of the atom
411 
412   real(dp), allocatable :: coretau(:)
413    ! coretau(mesh_size)
414    ! Gives the kinetic energy density of the atom
415 
416   real(dp), allocatable :: dij0(:)
417    ! dij0(lmn2_size)
418    ! Part of the Dij term (non-local operator) completely
419    ! calculated in the atomic data part
420 
421   real(dp), allocatable :: dltij(:)
422    ! dltij(lmn2_size)
423    ! Factor used to compute sums over klmn=(ilmn,jlmn)
424    ! ((ilmn,ilmn) term has to be added once)
425    ! dltij(klmn)=1 if ilmn=jlmn, else dltij(klmn)=2
426 
427   real(dp), allocatable :: dshpfunc(:,:,:)
428    ! shapefunc(mesh_size,l_size,4)
429    ! Gives the 4 first derivatives of  radial shape function
430    ! for each l component; used only if shape_type=-1
431 
432   real(dp), allocatable :: eijkl(:,:)
433    ! eijkl(lmn2_size,lmn2_size)
434    ! Hartree kernel for the on-site terms (E_hartree=Sum_ijkl[rho_ij rho_kl e_ijkl])
435    ! Used for Hartree and/or Fock contributions
436 
437   real(dp), allocatable :: eijkl_sr(:,:)
438    ! eijkl_sr(lmn2_size,lmn2_size)
439    ! Screened Hartree kernel for the on-site terms (E_hartree=Sum_ijkl[rho_ij rho_kl e_ijkl_sr])
440    ! Used for screened Fock contributions
441 
442   real(dp), allocatable :: euijkl(:,:,:,:,:)
443    ! euijkl(3,lmn_size,lmn_size,lmn_size,lmn_size)
444    ! PAW+U kernel for the on-site terms ( E_PAW+U = 0.5 * Sum_ijkl Sum_s1s2 [rho_ij^s1 rho_kl^s2 euijkl(s1,s2)] )
445    ! Contrary to eijkl and eijkl_sr, euijkl is not invariant with respect to the permutations i <--> j or k <--> l
446    ! However, it is still invariant with respect to the permutation i,k <--> j,l, see pawpuxinit.F90
447    ! Also, it depends on two spin indexes
448    ! Used for PAW+U contributions
449 
450   real(dp), allocatable :: euij_fll(:)
451    ! euij_fll(lmn2_size)
452    ! Double counting part of the PAW+U kernel in the "fully localized limit".This term is only linear with respect to rho_ij,
453    ! while euijkl is quadratic.
454    ! Used for PAW+U contributions
455 
456   real(dp), allocatable :: ex_cvij(:)
457   ! ex_cvij(lmn2_size))
458   ! Onsite exact_exchange matrix elements for core-valence interactions of the Fock operator
459 
460   real(dp), allocatable :: fk(:,:)
461    ! fk(6,4)
462    ! Slater integrals used for local exact exchange
463 
464   real(dp), allocatable :: gammaij(:)
465    ! gammaij(lmn2_size)
466    ! background contribution from the densities
467 
468   real(dp), allocatable :: gnorm(:)
469    ! gnorm(l_size)
470    ! Give the the normalization factor of each radial shape function
471 
472   real(dp), allocatable :: kij(:)
473   ! kij(lmn2_size))
474   ! Onsite matrix elements <phi|\kinetic|phj>-<tphi|\kinetic|tphj>
475 
476   real(dp), allocatable :: nabla_ij(:,:,:)
477    ! nabla_ij(3,lmn_size,lmn_size)
478    ! Onsite matrix elements <phi|\nabla|phj>-<tphi|\nabla|tphj>
479 
480   real(dp), allocatable :: nabla_im_ij(:,:,:)
481    ! nabla_im_ij(3,lmn_size,lmn_size)
482    ! Imaginary part of onsite matrix elements <phi|\nabla|phj>-<tphi|\nabla|tphj>
483    ! Used in case of core spinor wave functions
484 
485   real(dp), allocatable :: nablaphi(:,:)
486    ! nablaphi(partialwave_mesh_size, basis_size)
487    ! store the results of dphi/dr-(1/r)phi
488 
489   real(dp), allocatable :: phi(:,:)
490    ! phi(partialwave_mesh_size, basis_size)
491    ! Gives the paw electron wavefunctions on the radial grid
492 
493   real(dp), allocatable :: phiphj(:,:)
494    ! phiphj(mesh_size,ij_size)
495    ! Useful product Phi(:,i)*Phi(:,j)
496 
497   real(dp), allocatable :: phiphjint(:)
498    ! phiphjint(ij_proj)
499    ! Integration of Phi(:,i)*Phi(:,j) for DFT+U/local exact-exchange occupation matrix
500 
501   real(dp), allocatable :: ph0phiint(:)
502    ! ph0phjint(ij_proj)
503    ! Integration of Phi(:,1)*Phi(:,j) for LDA+DMFT projections
504 
505   real(dp), allocatable :: qgrid_shp(:)
506    ! qgrid_shp(mqgrid_shp)
507    ! Grid of points in reciprocal space on which the shape functions are given
508 
509   real(dp), allocatable :: qijl(:,:)
510    ! qijl(l_size**2,lmn2_size)
511    ! The qijl are the moments of the charge density difference between
512    ! the AE and PS partial wave for each channel (i,j). They take part
513    ! to the building of the compensation charge
514 
515   real(dp), allocatable :: rad_for_spline(:)
516    ! rad_for_spline(mesh_size)
517    ! Radial mesh used to spline quantities on radial mesh;
518    ! Allocated and used only when
519    !     shape_type=-1 (numerical shape function)
520    !  or usedvloc=1 (use of vloc derivative)
521 
522   real(dp), allocatable :: rhoij0(:)
523    ! rhoij0(lmn2_size)
524    ! Initial guess for rhoij
525 
526   real(dp), allocatable :: shape_alpha(:,:)
527    ! shape_alpha(2,l_size)
528    ! Alpha_i parameters in Bessel shapefunctions (shape_type=3)
529 
530   real(dp), allocatable :: shape_q(:,:)
531    ! shape_q(2,l_size)
532    ! Q_i parameters in Bessel shapefunctions (shape_type=3)
533 
534   real(dp), allocatable :: shapefunc(:,:)
535    ! shapefunc(mesh_size,l_size)
536    ! Gives the normalized radial shape function for each l component
537 
538   real(dp), allocatable :: shapefncg(:,:,:)
539    ! shapefncg(mqgrid_shp,2,l_size)
540    ! Gives the spherical Fourier transform of the radial shape function
541    ! for each l component (for each qgrid_shp(i)) + second derivative
542 
543   real(dp), allocatable :: sij(:)
544    ! sij(lmn2_size)
545    ! Nonlocal part of the overlap operator
546 
547   real(dp), allocatable :: tcoredens(:,:)
548    ! tcoredens(core_mesh_size,1)
549    ! Gives the pseudo core density of the atom
550    ! In PAW+WVL:
551    !  tcoredens(core_mesh_size,2:6)
552    !  are the first to the fifth derivatives of the pseudo core density.
553 
554  real(dp), allocatable :: tcoretau(:)
555    ! tcoretau(coretau_mesh_size)
556    ! Gives the pseudo core kinetic energy density of the atom
557 
558   real(dp), allocatable :: tcorespl(:,:)
559    ! tcorespl(mqgrid,2)
560    ! Gives the pseudo core density in reciprocal space on a regular grid
561 
562   real(dp), allocatable :: tcoretauspl(:,:)
563    ! tcoretauspl(mqgrid,2)
564    ! Gives the pseudo kinetic core density in reciprocal space on a regular grid
565 
566   real(dp), allocatable :: tnablaphi(:,:)
567    ! tphi(partialwave_mesh_size,basis_size)
568    ! Gives, on the radial grid, the paw atomic pseudowavefunctions
569 
570   real(dp), allocatable :: tphi(:,:)
571    ! tphi(partialwave_mesh_size,basis_size)
572    ! Gives, on the radial grid, the paw atomic pseudowavefunctions
573 
574   real(dp), allocatable :: tphitphj(:,:)
575    ! tphitphj(mesh_size,ij_size)
576    ! Useful product tPhi(:,i)*tPhi(:,j)
577 
578   real(dp), allocatable :: tproj(:,:)
579    ! non-local projectors
580 
581   real(dp), allocatable :: tvalespl(:,:)
582    ! tvalespl(mqgrid,2)
583    ! Gives the pseudo valence density in reciprocal space on a regular grid
584 
585   real(dp), allocatable :: Vee(:,:,:,:)
586    ! PAW+U:
587    ! Screened interaction matrix deduced from U and J parameters
588    ! computed on the basis of orbitals on which U acts.
589 
590   real(dp), allocatable :: Vex(:,:,:,:,:)
591    ! Local exact-exchange:
592    ! Screened interaction matrix deduced from calculation of Slater integrals
593    ! computed on the basis of orbitals on which local exact exchange acts.
594 
595   real(dp), allocatable :: vhtnzc(:)
596    ! vhtnzc(mesh_size)
597    ! Hartree potential for pseudized Zc density, v_H[\tilde{n}_{Zc}]
598    ! read in from PAW file
599 
600   real(dp), allocatable :: VHnZC(:)
601    ! VHnZC(mesh_size)
602    ! Hartree potential for Zc density, v_H[n_{Zc}]
603    ! constructed from core density in PAW file (see psp7in.F90)
604 
605   real(dp), allocatable :: vminushalf(:)
606    ! vminushalf(mesh_size)
607    ! External potential for LDA minus half calculation
608    ! read in from PAW file
609 
610   real(dp), allocatable :: zioneff(:)
611    ! zioneff(ij_proj)
612    ! "Effective charge"*n "seen" at r_paw, deduced from Phi at r_paw, n:
613    ! pricipal quantum number
614    ! good approximation to model wave function outside PAW-sphere through
615 
616  end type pawtab_type
617 
618  public :: pawtab_free         ! Free memory
619  public :: pawtab_nullify      ! Nullify content
620  public :: pawtab_get_lsize    ! Get the max. l for a product of 2 partial waves
621  public :: pawtab_set_flags    ! Set the value of the internal flags
622  public :: pawtab_print        ! Printout of the object.
623  public :: pawtab_bcast        ! MPI broadcast the object
624 
625  interface pawtab_nullify
626    module procedure pawtab_nullify_0D
627    module procedure pawtab_nullify_1D
628  end interface pawtab_nullify
629 
630  interface pawtab_free
631    module procedure pawtab_free_0D
632    module procedure pawtab_free_1D
633  end interface pawtab_free
634 
635  interface pawtab_set_flags
636    module procedure pawtab_set_flags_0D
637    module procedure pawtab_set_flags_1D
638  end interface pawtab_set_flags

m_pawtab/wvlpaw_allocate [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  wvlpaw_allocate

FUNCTION

  Allocate (if necessary) and nullify content of a wvlpaw pointer

SIDE EFFECTS

  wvlpaw<type(wvlpaw_type)>=datastructure to be allocated.

SOURCE

2946 subroutine wvlpaw_allocate(wvlpaw)
2947 
2948 !Arguments ------------------------------------
2949  type(wvlpaw_type),pointer :: wvlpaw
2950 
2951 ! *************************************************************************
2952 
2953  !@wvlpaw_type
2954 
2955  if (.not.associated(wvlpaw)) then
2956    LIBPAW_DATATYPE_ALLOCATE(wvlpaw,)
2957    call wvlpaw_nullify(wvlpaw)
2958  end if
2959 
2960  wvlpaw%npspcode_init_guess=10
2961 
2962 end subroutine wvlpaw_allocate

m_pawtab/wvlpaw_free [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  wvlpaw_free

FUNCTION

  Deallocate arrays and nullify flags in a wvlpaw structure

SIDE EFFECTS

  wvlpaw<type(wvlpaw_type)>=datastructure to be destroyed.
  All allocated arrays are deallocated.

SOURCE

2980 subroutine wvlpaw_free(wvlpaw)
2981 
2982 !Arguments ------------------------------------
2983  type(wvlpaw_type),pointer :: wvlpaw
2984 
2985 ! *************************************************************************
2986 
2987  !@wvlpaw_type
2988 
2989  if (.not.associated(wvlpaw)) return
2990 
2991  if(allocated(wvlpaw%pngau)) then
2992    LIBPAW_DEALLOCATE(wvlpaw%pngau)
2993  end if
2994  if(allocated(wvlpaw%parg)) then
2995    LIBPAW_DEALLOCATE(wvlpaw%parg)
2996  end if
2997  if(allocated(wvlpaw%pfac)) then
2998    LIBPAW_DEALLOCATE(wvlpaw%pfac)
2999  end if
3000 
3001  wvlpaw%npspcode_init_guess=0
3002  wvlpaw%ptotgau=0
3003 
3004  call wvlpaw_rholoc_free(wvlpaw%rholoc)
3005 
3006  LIBPAW_DATATYPE_DEALLOCATE(wvlpaw)
3007 
3008 end subroutine wvlpaw_free

m_pawtab/wvlpaw_nullify [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  wvlpaw_nullify

FUNCTION

  Nullify flags in a wvlpaw structure

SIDE EFFECTS

  wvlpaw=datastructure to be nullified

SOURCE

3025 subroutine wvlpaw_nullify(wvlpaw)
3026 
3027 !Arguments ------------------------------------
3028  type(wvlpaw_type),pointer :: wvlpaw
3029 
3030 ! *************************************************************************
3031 
3032  !@wvlpaw_type
3033  if (.not.associated(wvlpaw)) return
3034 
3035  wvlpaw%npspcode_init_guess=0
3036  wvlpaw%ptotgau=0
3037 
3038  call wvlpaw_rholoc_nullify(wvlpaw%rholoc)
3039 
3040 end subroutine wvlpaw_nullify

m_pawtab/wvlpaw_rholoc_free [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  wvlpaw_rholoc_free

FUNCTION

  Deallocate arrays and nullify flags in a wvlpaw%rholoc structure

SIDE EFFECTS

  wvlpaw_rholoc<type(wvlpaw_rholoc_type)>=datastructure to be destroyed.
  All allocated arrays are deallocated.

SOURCE

3058 subroutine wvlpaw_rholoc_free(wvlpaw_rholoc)
3059 
3060 !Arguments ------------------------------------
3061  type(wvlpaw_rholoc_type),intent(inout) :: wvlpaw_rholoc
3062 
3063 ! *************************************************************************
3064 
3065  !@wvlpaw_rholoc_type
3066 
3067  if(allocated(wvlpaw_rholoc%d)) then
3068    LIBPAW_DEALLOCATE(wvlpaw_rholoc%d)
3069  end if
3070  if(allocated(wvlpaw_rholoc%rad)) then
3071    LIBPAW_DEALLOCATE(wvlpaw_rholoc%rad)
3072  end if
3073 
3074  wvlpaw_rholoc%msz=0
3075 
3076 end subroutine wvlpaw_rholoc_free

m_pawtab/wvlpaw_rholoc_nullify [ Functions ]

[ Top ] [ m_pawtab ] [ Functions ]

NAME

  wvlpaw_rholoc_nullify

FUNCTION

  Nullify flags in a wvlpaw%rholoc structure

SIDE EFFECTS

  wvlpaw_rholoc<type(wvlpaw_rholoc_type)>=datastructure to be nullified.

SOURCE

3093 subroutine wvlpaw_rholoc_nullify(wvlpaw_rholoc)
3094 
3095 !Arguments ------------------------------------
3096  type(wvlpaw_rholoc_type),intent(inout) :: wvlpaw_rholoc
3097 
3098 ! *************************************************************************
3099 
3100  !@wvlpaw_rholoc_type
3101 
3102  wvlpaw_rholoc%msz=0
3103 
3104 end subroutine wvlpaw_rholoc_nullify

m_pawtab/wvlpaw_rholoc_type [ Types ]

[ Top ] [ m_pawtab ] [ Types ]

NAME

 wvlpaw_rholoc_type

FUNCTION

 Objects for WVL+PAW

SOURCE

47  type,public :: wvlpaw_rholoc_type
48 
49   integer :: msz
50 ! mesh size
51 
52   real(dp),allocatable :: d(:,:)
53 ! local rho and derivatives
54 
55   real(dp),allocatable :: rad(:)
56 ! radial mesh
57 
58  end type wvlpaw_rholoc_type
59 
60  public :: wvlpaw_rholoc_free    ! Free memory
61  public :: wvlpaw_rholoc_nullify ! Nullify content

m_pawtab/wvlpaw_type [ Types ]

[ Top ] [ m_pawtab ] [ Types ]

NAME

 wvlpaw_type

FUNCTION

 Objects for WVL+PAW

SOURCE

 75  type,public :: wvlpaw_type
 76 
 77 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
 78 ! declared in another part of ABINIT, that might need to take into account your modification.
 79 
 80 !Integer scalars
 81 
 82   integer :: npspcode_init_guess
 83    ! This is for the PAW-WVL case, only for the initial guess
 84 
 85   integer :: ptotgau
 86    ! total number of complex gaussians
 87    ! for tproj
 88 
 89   integer,allocatable :: pngau(:)
 90    ! number of complex gaussians per basis element
 91    ! for tproj
 92 
 93 !Real pointers
 94 
 95   real(dp),allocatable :: parg(:,:)
 96    !argument of Gaussians
 97 
 98   real(dp),allocatable :: pfac(:,:)
 99    !factors of Gaussians
100 
101 !Other scalars
102 
103   type(wvlpaw_rholoc_type) :: rholoc
104    ! local density
105    !   d(:,1): local rho
106    !   d(:,2): local rho 2nd-derivative
107    !   d(:,3): local pot
108    !   d(:,4): local pot 2nd-derivative
109 
110  end type wvlpaw_type
111 
112  public :: wvlpaw_allocate  ! Allocate memory
113  public :: wvlpaw_free   ! Free memory
114  public :: wvlpaw_nullify

m_pawtap/pawtab_get_lsize [ Functions ]

[ Top ] [ Functions ]

NAME

  pawtab_get_lsize

FUNCTION

  From an array of pawtab datastructures, get, for each atom, the value
  of "l_size" parameter.
  l_size is the maximum value of l accessible by a product of 2 partial waves;
  it may be cut by dtset%pawlcutd parameter

INPUTS

   [mpi_atmtab(:)]=--optional-- indexes of the atoms treated by current proc
   natom= number of atoms (may be a local or absolute number of atoms)
   typat(:)= list of atom types

OUTPUT

   l_size_atm(natom)=output array of l_size values (for each atom)

NOTES

  This function returns an allocatable integer array which may be allocated
  on the fly.

SOURCE

1354 subroutine pawtab_get_lsize(Pawtab,l_size_atm,natom,typat, &
1355 &                           mpi_atmtab) ! Optional argument
1356 
1357 !Arguments ------------------------------------
1358 !scalars
1359  integer,intent(in) :: natom
1360 !arrays
1361  integer,intent(in) :: typat(:)
1362  integer,optional,intent(in) :: mpi_atmtab(:)
1363  integer,allocatable,intent(inout) :: l_size_atm(:)
1364  type(pawtab_type),intent(in) :: pawtab(:)
1365 
1366 !Local variables-------------------------------
1367  integer :: ia,ityp,natom_typat
1368  character(len=100) :: msg
1369 
1370 ! *************************************************************************
1371 
1372  !@pawtab_type
1373 
1374  natom_typat=count(typat>0)
1375  if (size(pawtab)<maxval(typat)) then
1376    msg='error on pawtab size!'
1377    LIBPAW_BUG(msg)
1378  end if
1379 
1380  if (.not.allocated(l_size_atm)) then
1381    LIBPAW_ALLOCATE(l_size_atm,(natom))
1382  else if (size(l_size_atm)/=natom) then
1383    LIBPAW_DEALLOCATE(l_size_atm)
1384    LIBPAW_ALLOCATE(l_size_atm,(natom))
1385  end if
1386 
1387  if (natom==0) return
1388 
1389  if (natom==natom_typat) then
1390 
1391 !First case: sequential mode
1392    do ia=1,natom
1393      ityp=typat(ia)
1394      l_size_atm(ia)=pawtab(ityp)%lcut_size
1395    end do
1396 
1397  else
1398 
1399 !2nd case: parallel mode
1400    if (.not.present(mpi_atmtab)) then
1401      msg='optional args error!'
1402      LIBPAW_BUG(msg)
1403    end if
1404    do ia=1,natom
1405      ityp=typat(mpi_atmtab(ia))
1406      l_size_atm(ia)=pawtab(ityp)%lcut_size
1407    end do
1408 
1409  end if
1410 
1411 end subroutine pawtab_get_lsize