User Tools

Site Tools


atomicdata:specs:paw-xml

Differences

This shows you the differences between two versions of the page.

Link to this comparison view

Both sides previous revisionPrevious revision
Next revision
Previous revision
atomicdata:specs:paw-xml [2014/11/01 21:53] – Add PAWdataset example file (N.lda) Marc Torrentatomicdata:specs:paw-xml [2015/05/18 00:25] (current) – Change def of N_c (factor 2 was missing) Marc Torrent
Line 1: Line 1:
 +====== XML specification for atomic PAW datasets ======
 +
 +<HTML><br></HTML>
 +//Current version of the specification: **0.7**// \\ 
 +//Codes implementing the PAW XML format:
 +[[https://wiki.fysik.dtu.dk/gpaw|gpaw]],
 +[[http://www.abinit.org|abinit]],
 +[[http://users.wfu.edu/natalie/papers/pwpaw/man.html|atompaw]]//
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Introduction =====
 +<HTML><br></HTML>
 +
 +This page contains information about the PAW-XML data format for the atomic datasets
 +necessary for doing Projector Augmented-Wave calculations [(BLOCHL1994)].
 +We use the term //dataset// instead of //pseudo-potential// because
 +the PAW method is not a pseudo-potential method.
 +
 +An example XML file for nitrogen PAW dataset using LDA can be seen here: {{atomicdata:specs:n.lda.xml| n.lda.xml}}.
 +
 +<HTML><blockquote>
 +<b>Note</b>: <u>Hartree</u> atomic units are used in the XML file (&hbar; = <i>m</i> = <i>e</i> = 1).
 +</blockquote></HTML>
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== What defines a dataset? =====
 +<HTML><br></HTML>
 +
 +The following quantities define a minimum PAW dataset (the notation from Ref. [(BLOCHL2003)] is used here):
 +
 +^Quantity                    ^Description                           ^
 +|Z                           |Atomic number                         |
 +|$E_\text{XC}[n]$            |Exchange-correlation functional       |
 +|$E^\text{kin}_c$            |Kinetic energy of the core electrons  |
 +|$g_{\ell m}(\mathbf{r})$    |Shape function for compensation charge|
 +|$n_c(r)$                    |All-electron core density             |
 +|$\tilde{n}_c(r)$            |Pseudo electron core density          |
 +|$\tilde{n}_v(r)$            |Pseudo electron valence density       |
 +|$\bar{v}(r)$                |Zero potential                        |
 +|$\phi_i(\mathbf{r})$        |All-electron partial waves            |
 +|$\tilde{\phi}_i(\mathbf{r})$|Pseudo partial waves                  |
 +|$\tilde{p}_i(\mathbf{r})$   |Projector functions                   |
 +|$\Delta E^\text{kin}_{ij}$  |Kinetic energy differences            |
 +
 +\\
 +The following quantities can be optionally provided:
 +
 +^Quantity                    ^Description                                                   ^
 +|$r_{PAW}$                   |Radius of the PAW augmentation region (max. of matching radii)|
 +|$v_H[\tilde{n}_{Zc}](r)$    |Kresse-Joubert local ionic pseudopotential                    |
 +|$\tilde{Q}_{ij}(\mathbf{r})$|State-dependent shape function for compensation charge        |
 +|$\tau_c(r)$                 |Core kinetic energy density                                   |
 +|$\tilde{\tau}_c(r)$         |Pseudo core kinetic energy density                            |
 +|$X^{\text{core-core}}$      |Core-core contribution to exact exchange                      |
 +|$X_{ij}^{\text{core-val}}$  |Core-valence exact-exchange correction matrix                 |
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Specification of the elements =====
 +<HTML><br></HTML>
 +
 +An element looks like this:
 +<code xml>
 +<name> ... </name>
 +</code>
 +
 +or for an empty element:
 +<code xml>
 +<name/>
 +</code>
 +
 +<HTML><blockquote>
 +<b>Tip</b>: An XML-tutorial can be found
 +<a href="http://www.w3schools.com/xml/default.asp">here</a>.
 +</blockquote></HTML>
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== The header =====
 +<HTML><br></HTML>
 +
 +The first two lines should look like this:
 +<code xml>
 +<?xml version="1.0"?>
 +<paw_dataset version="0.7">
 +</code>
 +
 +The first line must be present in all XML files.
 +Everything else is put inside an element with name ''%%paw_dataset%%'', and this element has an attribute called ''%%version%%''.
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== A comment =====
 +<HTML><br></HTML>
 +
 +It is recommended to put a comment giving the units and a link to this web page:
 +<code xml>
 +<!-- Nitrogen dataset for the Projector Augmented Wave method. -->
 +<!-- Units: Hartree and Bohr radii.                            -->
 +<!-- http://www.where.org/paw_dataset.html                     -->
 +</code>
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== The ''atom'' element =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +<atom symbol="N" Z="7" core="2" valence="5"/>
 +</code>
 +
 +The ''%%atom%%'' element has attributes ''%%symbol%%'', ''%%Z%%'', ''%%core%%'' and ''%%valence%%'' (chemical symbol, atomic number, number of core electrons and number of valence electrons).
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Exchange-correlation =====
 +<HTML><br></HTML>
 +
 +The ''%%xc_functional%%'' element defines the exchange-correlation functional used for generating the dataset. It has the two attributes ''%%type%%'' and ''%%name%%''.
 +
 +The ''%%type%%'' attribute can be ''%%LDA%%'', ''%%GGA%%'', ''%%MGGA%%'' or ''%%HYB%%''.
 +
 +The ''%%name%%'' attribute designates the exchange-correlation functional and can be specified in the following ways:
 +
 +<HTML><li></HTML>
 +Taking the names from the [[http://www.tddft.org/programs/octopus/wiki/index.php/Libxc:manual#Available_functionals|LibXC]] library. \\ The correlation and exchange names are stripped from their ''%%XC_%%'' part and combined with a ''%%+%%''-sign. \\ Here is an example for an LDA functional:<code xml><xc_functional type="LDA", name="LDA_X+LDA_C_PW"/></code> and this is what PBE will look like:<code xml><xc_functional type="GGA", name="GGA_X_PBE+GGA_C_PBE"/></code>
 +<HTML></li></HTML>
 +
 +<HTML><li></HTML>
 +Using one of the following pre-defined aliases:
 +
 +^''%%type%%''^''%%name%%''  ^LibXC equivalent                   ^Reference                                                                      ^
 +|''%%LDA%%'' |''%%PW%%''    |''%%LDA_X+LDA_C_PW%%''             |LDA exchange; Perdew, Wang, //PRB// **45**, 13244 (1992)                       |
 +|''%%GGA%%'' |''%%PW91%%''  |''%%GGA_X_PW91+GGA_C_PW91%%''      |Perdew et al //PRB// **46**, 6671 (1992)                                       |
 +|''%%GGA%%'' |''%%PBE%%''   |''%%GGA_X_PBE+GGA_C_PBE%%''        |Perdew, Burke, Ernzerhof, //PRL// **77**, 3865 (1996)                          |
 +|''%%GGA%%'' |''%%RPBE%%''  |''%%GGA_X_RPBE+GGA_C_PBE%%''       |Hammer, Hansen, Nørskov, //PRB// **59**, 7413 (1999)                           |
 +|''%%GGA%%'' |''%%revPBE%%''|''%%GGA_X_PBE_R+GGA_C_PBE%%''      |Zhang, Yang, //PRL// **80**, 890 (1998)                                        |
 +|''%%GGA%%'' |''%%PBEsol%%''|''%%GGA_X_PBE_SOL+GGA_C_PBE_SOL%%''|Perdew et al, //PRL// **100**, 136406 (2008)                                   |
 +|''%%GGA%%'' |''%%AM05%%''  |''%%GGA_X_AM05+GGA_C_AM05%%''      |Armiento, Mattsson, //PRB// **72**, 085108 (2005)                              |
 +|''%%GGA%%'' |''%%BLYP%%''  |''%%GGA_X_B88+GGA_C_LYP%%''        |Becke, //PRA// **38**, 3098 (1988); Lee, Yang, Parr, //PRB// **37**, 785 (1988)|
 +
 +Examples:
 +<code xml>
 +<xc_functional type="LDA", name="PW"/>
 +</code>
 +<code xml>
 +<xc_functional type="GGA", name="PBE"/>
 +</code>
 +<HTML></li></HTML>
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Generator =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +<generator type="scalar-relativistic" name="MyGenerator-2.0">
 + Frozen core: [He]
 +</generator>
 +</code>
 +
 +This element contains //character data// describing in words how the dataset was generated.
 +The ''%%type%%'' attribute must be one of: ''%%non-relativistic%%'', ''%%scalar-relativistic%%'' or ''%%relativistic%%''.
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Energies =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +<ae_energy kinetic="53.777460" xc="-6.127751"
 +           electrostatic="-101.690410" total="-54.040701"/>
 +<core_energy kinetic="43.529213"/>
 +</code>
 +The kinetic energy of the core electrons, $E^\text{kin}_c$, is used in the PAW method.
 +The other energies are convenient to have for testing purposes and can also be useful for checking
 +the quality of the underlying atomic calculation.
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Valence states =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +<valence_states>
 +  <state n="2" l="0" f="2"  rc="1.10" e="-0.6766" id="N-2s"/>
 +  <state n="2" l="1" f="3"  rc="1.10" e="-0.2660" id="N-2p"/>
 +  <state       l="0"        rc="1.10" e=" 0.3234" id="N-s1"/>
 +  <state       l="1"        rc="1.10" e=" 0.7340" id="N-p1"/>
 +  <state       l="2"        rc="1.10" e=" 0.0000" id="N-d1"/>
 +</valence_states>
 +</code>
 +
 +The ''%%valence_states%%'' element contains several ''%%state%%'' elements,
 +defined by a unique ''%%id%%'' as well as ''%%l%%'' and ''%%n%%'' quantum numbers.
 +For each of them it is also required to provide the energy ''%%e%%'', the occupation ''%%f%%''
 +and the matching radius of the partial waves ''%%rc%%''.\\ 
 +The number of ''%%state%%'' elements determines the size of the partial wave basis. It is equal to the number of [[paw-xml#Radial functions|radial functions]] (radial parts of the $\phi_i$,
 +$\tilde{\phi}_i$ and $\tilde{p}_i$) and is noted $n_{waves}$ in the rest of this document.
 +
 +For this dataset, the first two lines describe bound eigenstates with occupation numbers
 +and principal quantum numbers. Notice, that the three additional unbound states should
 +have no ''%%f%%'' and ''%%n%%'' attributes.
 +In this way, we know that only the first two bound states (with ''%%f%%'' and ''%%n%%'' attributes)
 +should be used for constructing an initial guess for the wave functions.
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Radial grids =====
 +<HTML><br></HTML>
 +
 +There can be one or more definitions of radial grids.
 +
 +Example:
 +<code xml>
 +<radial_grid eq="r=d*i" d="0.1" istart="0" iend="9" id="g1">
 +  <values>
 +    0.0 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9
 +  </values>
 +  <derivatives>
 +    0.1 0.1 0.1 0.1 0.1 0.1 0.1 0.1 0.1 0.1
 +  </derivatives>
 +</radial_grid>
 +</code>
 +
 +This defines one radial grid as $r_i = di$ where $i$ runs from 0 to 9.
 +Inside the ''%%<radial_grid>%%'' element we have the 10 values of $r_i$ followed by the 10 values
 +of the derivatives $dr_i/di$.
 +
 +All functions (densities, potentials, ...) that use this grid are given as 10 numbers defining
 +the radial part of the function. The radial part of the function must be multiplied by a spherical
 +harmonics: $f_{\ell m}(\mathbf{r}) = f_\ell(r) Y_{\ell m}(\theta, \phi)$.
 +
 +Each radial grid has a unique id:
 +<code xml>
 +<radial_grid eq="r=d*i" d="0.01" istart="0" iend="99" id="lin">
 +<radial_grid eq="r=a*exp(d*i)" a="1.056e-4" d="0.05" istart="0" iend="249" id="log">
 +</code>
 +
 +and each numerical function must refer to one of these ids:
 +<code xml>
 +<function grid="lin">
 +  ... ...
 +</function>
 +</code>
 +
 +In this example, the ''%%function%%'' element should contain 100 numbers ($i$ = 0, ..., 99).
 +Each number must be separated by a ''%%<newline>%%'' character or by one or more ''%%<tab>%%'''s
 +or ''%%<space>%%'''s (no commas).
 +For numbers with scientific notation, use this format: ''%%1.23456e-5%%'' or ''%%1.23456E-5%%''
 +and not ''%%1.23456D-5%%''.
 +
 +A program can read the values for $r_i$ and $dr_i/di$ from the file or evaluate them from the ''%%eq%%''
 +and associated parameter attributes. There are currently six types of radial grids:
 +
 +^''%%eq%%''               ^parameters             ^
 +|''%%r=d*i%%''            |''%%d%%''              |
 +|''%%r=a*exp(d*i)%%''     |''%%a%%'' and ''%%d%%''|
 +|''%%r=a*(exp(d*i)-1)%%'' |''%%a%%'' and ''%%d%%''|
 +|''%%r=a*i/(1-b*i)%%''    |''%%a%%'' and ''%%b%%''|
 +|''%%r=a*i/(n-i)%%''      |''%%a%%'' and ''%%n%%''|
 +|''%%r=(i/n+a)^5/a-a^4%%''|''%%a%%'' and ''%%n%%''|
 +
 +The ''%%istart%%'' and ''%%iend%%'' attributes indicating the range of $i$ should always be present.
 +
 +Although it is possible to define as radial grids as desired,
 +it is recommended to minimize the number of grids in the dataset.
 +
 +\\ \\  
 +<HTML><hr></HTML>
 +===== Shape function for the compensation charge =====
 +<HTML><br></HTML>
 +
 +The general formulation of the compensation charge uses an expansion over the partial waves //ij// and the spherical harmonics:
 +
 +$$\sum_{\ell m} C_{\ell m \ell_i m_i \ell_j m_j} \hat{Q}^{\ell}_{i j}(r) Y_{\ell m}(\theta, \phi),$$
 +
 +where $C_{\ell m \ell_i m_i \ell_j m_j}$ is a //Gaunt coefficient//.
 +
 +The standard expression [(BLOCHL1994)] for the //shape function// $\hat{Q}^{\ell}_{i j}(\mathbf{r})$
 +is a product of the multipole moment $Q^{\ell}_{i j}$ and a shape function $g_l(r)$:
 +
 +$$\hat{Q}^{\ell}_{i j}(r) = Q^{\ell}_{i j} g_\ell(r),$$
 +
 +
 +Several formulations [(BLOCHL1994)] [(HOLZWARTH2001)] define $g_\ell(r) \propto r^\ell k(r)$,
 +where $k(r)$ is an $l$-independent shape function:
 +
 +^''%%type%%'' ^parameters                 ^$k(r)$                           ^
 +|''%%gauss%%''|''%%rc%%''                 |$\exp(-(r/r_c)^2)$               |
 +|''%%sinc%%'' |''%%rc%%''                 |$[\sin(\pi r/r_c)/(\pi r/r_c)]^2$|
 +|''%%exp%%''  |''%%rc%%'' and ''%%lamb%%''|$\exp(-(r/r_c)^\lambda) $        |
 +
 +Example:
 +<code xml>
 +<shape_function type="gauss" rc="3.478505426185e-01">
 +</code>
 +
 +Another formulation [(KRESSE1999)] defines directly $g_l(r)$:
 +
 +^''%%type%%''  ^parameters       ^$g_l(r)$                                       ^
 +|''%%bessel%%''|''%%rc%%''       |$\sum_{i=1}^2 \alpha_i^\ell j_\ell(q_i^\ell r)$|
 +
 +For ''%%bessel%%'' the four parameters ($\alpha_1^l$, $q_1^l$, $\alpha_2^l$ and $q_2^l$)
 +must be determined from ''%%rc%%'' for each value of $l$ as described in [(KRESSE1999)].
 +
 +Example:
 +<code xml>
 +<shape_function type="bessel" rc="3.478505426185e-01">
 +</code>
 +
 +There is also a more general formulation where $\hat{Q}^{\ell}_{i j}(r)$ is given in a numerical form.
 +Several //shape functions// can be set (with the ''%%<shape_function>%%'' tag),
 +depending on $l$ and/or combinations of partial waves (specified using the optional ''%%state1%%''
 +and ''%%state2%%'' attributes). See for instance section II.C of [(LAASONEN1993)].
 +
 +Example 1, defining numerically $g_l(r)$ in $\hat{Q}^{\ell}_{i j}(r)=Q^{\ell}_{i j} g_\ell(r)$:
 +<code xml>
 +<shape_function type="numeric" l=0 grid="g1">
 + ... ... ...
 +</shape_function>
 +</code>
 +
 +Example 2, defining directly $\hat{Q}^{\ell}_{i j}(r)$ for states $i$=''%%N-2s%%'' and $j$=''%%N-2p%%'',
 +and $l$=0:
 +<code xml>
 +<shape_function type="numeric" l=0 state1="N-2s" state2="N-2p" grid="g1">
 + ... ... ...
 +</shape_function>
 +</code>
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Radial functions =====
 +<HTML><br></HTML>
 +
 +Continuing, we have now reached the //all-electron// (resp. //pseudo core//, //pseudo valence//) density:
 +<code xml>
 +<ae_core_density grid="g1">
 +  6.801207147443e+02 6.801207147443e+02 6.665042896724e+02
 +  ... ...
 +</ae_core_density>
 +
 +<pseudo_core_density rc="1.1" grid="g1">
 +  ... ...
 +</pseudo_core_density>
 +
 +<pseudo_valence_density rc="1.1" grid="g1">
 +  ... ...
 +</pseudo_valence_density>
 +</code>
 +
 +The numbers inside the ''%%ae_core_density%%'' (resp. ''%%pseudo_core_density%%'',
 +''%%pseudo_valence_density%%'') element defines the radial part of $n_c(\mathbf{r})$
 +(resp. $\tilde{n}_c(\mathbf{r})$, $\tilde{n}_v(\mathbf{r})$).
 +The radial part must be multiplied by $Y_{00} = (4\pi)^{-1/2}$ to get the full density
 +($Y_{00}n_c(\mathbf{r})$ should integrate to the number of core electrons).
 +The //pseudo core density// and the //pseudo valence// density are defined similarly and
 +also have a ''%%rc%%'' attribute specifying the matching radius.
 +
 +The ''%%ae_partial_wave%%'', ''%%pseudo_partial_wave%%'' and ''%%projector_function%%''
 +elements contain the radial parts of the $\phi_i(\mathbf{r})$,
 +$\tilde{\phi}_i(\mathbf{r})$ and $\tilde{p}_i(\mathbf{r})$
 +functions for the states listed in the ''%%valence_states%%'' element above
 +(five states in the nitrogen example).
 +All functions must have an attribute ''%%state="..."%%'' referring to one of the
 +states listed in the ''%%valence_states%%'' element:
 +<code xml>
 +<ae_partial_wave state="N-2s" grid="g1">
 +  -8.178800366898029e+00 -8.178246914143839e+00 -8.177654917302689e+00
 +  ... ...
 +</ae_partial_wave>
 +
 +<pseudo_partial_wave state="N-2s" grid="g1">
 +  ... ...
 +</pseudo_partial_wave>
 +
 +<projector_function state="N-2s" grid="g1">
 +  ... ...
 +</projector_function>
 +
 +<ae_partial_wave state="N-2p" grid="g1">
 +  ... ...
 +</ae_partial_wave>
 +</code>
 +
 +Remember that the radial part of these functions must be multiplied by a spherical harmonics:
 +$$\phi_i(\mathbf{r}) = \phi_i(r) Y_{\ell_i m_i}(\theta, \phi)$$
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Zero potential =====
 +<HTML><br></HTML>
 +
 +The zero potential, $\bar{v}$ (see section VI.D of [(BLOCHL1994)])
 +is defined similarly to the densities; the radial part must be multiplied by
 +$Y_{00} = (4\pi)^{-1/2}$ to get the full potential.
 +The ''%%zero_potential%%'' element has a ''%%rc%%'' attribute specifying the cut-off
 +radius of $\bar{v}(\mathbf{r})$:
 +<code xml>
 +<zero_potential rc="1.1" grid="g1">
 +  ... ...
 +</zero_potential>
 +</code>
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== The Kresse-Joubert formulation =====
 +<HTML><br></HTML>
 +
 +The Kresse-Joubert formulation of the PAW method [(KRESSE1999)]
 +is very similar to the original formulation of Blöchl [(BLOCHL1994)]
 +. However, the Kresse-Joubert formulation does not use $\bar{v}$ directly,
 +but indirectly through the local ionic pseudopotential, $v_H[\tilde{n}_{Zc}]$.
 +Therefore, the following transformation is necessary:
 +
 +$$v_H[\tilde{n}_{Zc}] = v_H[\tilde{n}_c +
 +(N_c - Z - \tilde{N}_c) g_{00} Y_{00}] + \bar{v} +
 +v_{xc}[\tilde{n}_v + \tilde{n}_c] -
 +v_{xc}[\tilde{n}_v + \tilde{n}_c +
 +       (N_v - \tilde{N}_v - \tilde{N}_c) g_{00} Y_{00}]$$
 +
 +where $N_c$ is the number of core electrons, $N_v$ is the number of valence electrons,
 +$\tilde{N}_c$ is the number of electrons contained in the pseudo core density and
 +$\tilde{N}_v$ is the number of electrons contained in the pseudo valence density.
 +The Hartree potential from the density n is defined as:
 +$$v_H[n](r_1) = 4\pi \int_0^\infty r_2^2 dr_2 \frac{n(r_2)}{r_>},$$
 +
 +where $r_>$ is the larger of $r_1$ and $r_2$.
 +
 +<HTML><blockquote>
 +<b>Note</b>:
 +In the Kresse-Joubert formulation, the symbol $\tilde{n}$ is used for
 +what we here call $\tilde{n}_v$ and in the Blöchl formulation,
 +we have $\tilde{n} = \tilde{n}_c + \tilde{n}_v$.
 +</blockquote></HTML>
 +
 +It is also possible to add an element ''%%kresse_joubert_local_ionic_pseudopotential%%''
 +that contains the $v_H[\tilde{n}_{Zc}](r)$ function directly, so that no conversion is necessary:
 +<code xml>
 +<kresse_joubert_local_ionic_pseudopotential rc="1.3" grid="log">
 +  ... ...
 +</kresse_joubert_local_ionic_pseudopotential>
 +</code>
 +
 +The ''%%kresse_joubert_local_ionic_pseudopotential%%'' element has a ''%%rc%%'' attribute
 +specifying the matching radius. This matching radius corresponds to the maximum of all
 +the matching radii used in the formalism.
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Kinetic energy differences =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +<kinetic_energy_differences>
 +  1.744042161013e+00 0.000000000000e+00 2.730637956456e+00
 +  ... ...
 +<kinetic_energy_differences>
 +</code>
 +
 +This element contains the symmetric $\Delta E^\text{kin}_{ij}$ matrix:
 +
 +$$\Delta E^\text{kin}_{ij} = \langle \phi_i | \hat{T} | \phi_j \rangle 
 +- \langle \tilde{\phi}_i | \hat{T} | \tilde{\phi}_j \rangle$$
 +
 +where $\hat{T}$ is the kinetic energy operator used by the generator.\\ 
 +With $n_{waves}$ partial waves (see [[paw-xml#Valence states|$n_{waves}$ definition]]), we have a $n_{waves} \times n_{waves}$ matrix listed as $n_{waves}^2$ numbers.
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Meta-GGA =====
 +<HTML><br></HTML>
 +
 +Datasets for use with MGGA functionals must also include information on the
 +//core kinetic energy density// and //pseudo core kinetic energy density// ;
 +the latters are defined with these two elements:
 +
 +<code xml>
 +<ae_core_kinetic_energy_density grid="g1"> 
 +  ... ...
 +</ae_core_kinetic_energy_density> 
 +<pseudo_core_kinetic_energy_density rc="1.1" grid="g1">
 +  ... ...
 +</pseudo_core_kinetic_energy_density>
 +</code>
 +
 +These densities are defined similarly to the core and valence densities (see above).
 +The ''%%pseudo_core_kinetic_energy_density%%'' element has a ''%%rc%%''
 +attribute specifying its matching radius.
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Exact exchange integrals =====
 +<HTML><br></HTML>
 +
 +The core-core contribution to the exact exchange energy $X^{\text{core-core}}$
 +and the symmetric core-valence PAW-correction matrix $X_{ij}^{\text{core-valence}}$
 +are given as:
 +
 +$$X^{\text{core-core}} = -\frac{1}{4}\sum_{cc'} \iint d\mathbf{r} d\mathbf{r}'
 +\frac{\phi_c(\mathbf{r})\phi_{c'}(\mathbf{r}) \phi_c(\mathbf{r}')\phi_{c'}(\mathbf{r}')}
 +{|\mathbf{r}-\mathbf{r}'|}$$
 +
 +$$X_{ij}^{\text{core-valence}} = -\frac{1}{2}\sum_c \iint d\mathbf{r} d\mathbf{r}'
 +\frac{\phi_i(\mathbf{r})\phi_c(\mathbf{r}) \phi_j(\mathbf{r}')\phi_c(\mathbf{r}')}
 +{|\mathbf{r}-\mathbf{r}'|}$$
 +
 +The $X_{ij}^{\text{core-valence}}$ coefficients depend only on pairs of the radial
 +basis functions $\phi_i(r)$ and can be evaluated by summing over radial
 +integrals times **3-j** symbols according to:
 +
 +$$X_{ij}^{\text{core-valence}} =
 +-\delta_{l_i l_j} \delta_{m_i m_j} \sum_{c L} \frac{N_c}{2}
 +{\begin{pmatrix}l_c & L & l_i \\ 0 & 0 & 0\end{pmatrix}}^2
 +\int r^2 dr \int {r'}^2 d{r'}
 +\frac{r^{L}_{<}}{r^{L+1}_{>}}
 +\phi_i(r) \phi_c(r) \phi_j(r') \phi_c(r')$$
 +
 +where\\ 
 +$N_{c}$ is the number of core electrons corresponding to $l_{c}$, namely $N_{c}=2(2l_{c}+1)$,\\ 
 +$r_>$ (resp. $r_<$) is the larger (resp. smaller) of $r$ and $r'$.
 +
 +$X^{\text{core-core}}$ can be specified in the ''%%core%%'' attribute of the ''%%<exact_exchange>%%''
 +element.\\ 
 +With $n_{waves}$ partial waves (see [[paw-xml#Valence states|$n_{waves}$ definition]]),
 +$X_{ij}^{\text{core-valence}}$ is a $n_{waves} \times n_{waves}$ matrix.
 +It can be specified as $n_{waves}^2$ numbers inside the ''%%<exact_exchange>%%'' element:
 +
 +<code xml>
 +<exact_exchange core="...">
 +  ... ...
 +</exact_exchange>
 +</code>
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Optional elements =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +<paw_radius rc="2.3456781234">
 +</code>
 +
 +Although not necessary, it may be helpful to provide the following item(s) in the dataset:
 +
 +  * Radius of the PAW augmentation region ''%%paw_radius%%'' \\ This radius defines the region (around the atom) outside which all pseudo quantities are equal to the all-electron ones. It is equal to the maximum of all the cut-off and matching radii. Note that -- for better lisibility -- the ''%%paw_radius%%'' element should be provided in the header of the file.
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== End of the dataset =====
 +<HTML><br></HTML>
 +
 +<code xml>
 +</paw_dataset>
 +</code>
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== How to use the datasets =====
 +<HTML><br></HTML>
 +
 +Most likely, the radial functions will be needed on some other type of radial grid than
 +the one used in the dataset. The idea is that one should read in the radial functions
 +and then transform them to the radial grids used by the specific implementation.
 +After the transformation, some sort of normalization may be necessary.
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== Plotting the radial functions =====
 +<HTML><br></HTML>
 +
 +The first 10-20 lines of the XML-datasets, should be pretty much human readable,
 +and should give an overview of what kind of dataset it is and how it was generated.
 +The remaining part of the file contains numerical data for all the radial functions.
 +To get an overview of these functions, you can extract that data with the
 +''%%pawxml.py%%'' program and then pass it on to your favorite plotting tool.
 +
 +<HTML><blockquote>
 +<b>Note</b>: The <span style="font-family: courier;">pawxml.py</span> program
 +is very primitive and is only included in order to demonstrates how to parse XML using
 +SAX from a Python program. Parsing XML from Fortran or C code with SAX should be similar.
 +</blockquote></HTML>
 +
 +<file python pawxml.py>
 +#!/usr/bin/env python
 +from __future__ import print_function
 +import os
 +import sys
 +import xml.sax
 +from optparse import OptionParser
 +from math import exp, log
 +
 +class Reader(xml.sax.handler.ContentHandler):
 +    def __init__(self, list, name, state):
 +        self.list = list
 +        self.name = name
 +        self.state = state
 +        self.grid = {}
 +        xml.sax.handler.ContentHandler.__init__(self)
 +    def read(self, filename):
 +        if filename.endswith('.gz'):
 +            import gzip
 +            source = gzip.open(filename)
 +        else:
 +            source = open(filename)
 +        xml.sax.parse(source, self)
 +    def startElement(self, name, attrs):
 +        if name == 'state' and self.list:
 +            print(attrs['id'], file=sys.stderr)
 +        elif name == 'radial_grid':
 +            istart = int(attrs['istart'])
 +            iend = int(attrs['iend'])
 +            ii = range(istart, iend + 1)
 +            id = attrs['id']
 +            eq = attrs['eq']
 +            if eq == 'r=a*exp(d*i)':
 +                a = float(attrs['a'])
 +                d = float(attrs['d'])
 +                self.grid[id] = [a * exp(d * i) for i in ii]
 +            elif eq == 'r=a*i/(n-i)':
 +                a = float(attrs['a'])
 +                n = float(attrs['n'])
 +                self.grid[id] = [a * i / (n - i) for i in ii]
 +            elif eq == 'r=a*(exp(d*i)-1)':
 +                a = float(attrs['a'])
 +                d = float(attrs['d'])
 +                self.grid[id] = [a * (exp(d * i) - 1.0) for i in ii]
 +            elif eq == 'r=d*i':
 +                d = float(attrs['d'])
 +                self.grid[id] = [d * i for i in ii]
 +            elif eq == 'r=(i/n+a)^5/a-a^4':
 +                a = float(attrs['a'])
 +                n = float(attrs['n'])
 +                self.grid[id] = [(i / n + a)**5 / a - a**4 for i in ii]
 +            else:
 +                raise RuntimeError('Unknown grid type: ' + attrs['eq'])
 +        elif name == self.name and attrs.get('state', None) == self.state:
 +            self.r = self.grid[attrs['grid']]
 +            self.data = ''
 +        else:
 +            self.data = None
 +    def characters(self, data):
 +        if self.data is not None:
 +            self.data += data
 +    def endElement(self, name):
 +        if self.data is not None:
 +            for r, x in zip(self.r, self.data.split()):
 +                print(r, x)
 +                
 +op = OptionParser(usage='%prog [options] setup[.gz]',
 +                      version='%prog 0.2')
 +op.add_option('-x', '--extract',
 +                  help='Function to extract.',
 +                  metavar='<name>')
 +op.add_option('-s', '--state',
 +                  help='Select valence state.',
 +                  metavar='<channel>')
 +op.add_option('-l', '--list', action='store_true',
 +                  help='List valence states.')
 +options, args = op.parse_args()
 +if len(args) != 1:
 +        op.error('incorrect number of arguments')
 +
 +reader = Reader(options.list, options.extract, options.state)
 +reader.read(args[0])
 +</file>
 +
 +Usage:
 +
 +It works like this:
 +<code>
 +$ pawxml.py [options] dataset[.gz]
 +</code>
 +
 +Options:
 +|''%%--version%%''                     |Show program's version number and exit.|
 +|''%%-h, --help%%''                    |Show this help message and exit.       |
 +|''%%-x <name>, --extract=<name>%%''   |Function to extract.                   |
 +|''%%-s<channel>, --state=<channel>%%''|Select valence state.                  |
 +|''%%-l, --list%%''                    |List valence states                    |
 +
 +Examples:
 +<code>
 +[~]$ pawxml.py -x pseudo_core_density N.LDA | xmgrace -
 +[~]$ pawxml.py -x ae_partial_wave -s N2p N.LDA > N.ae.2p 
 +[~]$ pawxml.py -x pseudo_partial_wave -s N2p N.LDA > N.ps.2p 
 +[~]$ xmgrace N.??.2p
 +</code>
 +
 +\\ \\ 
 +<HTML><hr></HTML>
 +===== References =====
 +<HTML><br></HTML>
 +
 +[(BLOCHL1994>[[http://journals.aps.org/prb/abstract/10.1103/PhysRevB.50.17953 | P. E. Blöchl, Projector Augmented-Wave method, Phys. Rev. B 50, 17953-19979 (1994) ]] )]
 +[(BLOCHL2003>[[http://www.ias.ac.in/matersci/bmsjan2003/Bms3867.pdf | P. E. Blöchl, C. J. Forst and J. Schimpl, Projector Augmented-Wave method: Ab initio molecular dynamics with full wave functions, Bulletin of Materials Science 26, 33-41 (2003) ]] )]
 +[(KRESSE1999>[[http://journals.aps.org/prb/abstract/10.1103/PhysRevB.59.1758 | G. Kresse and D. Joubert, Form ultrasoft pseudopotentials to the projector augmented-wave method, Phys. Rev. B 59, 1758-1775 (1999) ]] )]
 +[(HOLZWARTH2001>[[http://www.sciencedirect.com/science/article/pii/S0010465500002447 | N. A. W. Holzwarth, A. R. Tackett, and G. E. Matthews, A Projector Augmented Wave (PAW) code for electronics structure calculations: Part I atompaw for generating atom-centered functions, Computer Physics Communications 135, 329-347 (2001) ]] )]
 +[(LAASONEN1993>[[http://journals.aps.org/prb/abstract/10.1103/PhysRevB.47.10142 | K. Laasonen, A. Pasquarello, R. Car, C. Lee and D. Vanderbilt, Car-Parrinello molecular dynamics with Vanderbilt ultrasoft pseudopotentials, Phys. Rev. B 47, 10142-10153 (1993) ]] )]
 +
 +<refnotes>
 +refnote-id           : 1
 +reference-format     : []
 +reference-base       : text
 +multi-ref-id         : note
 +note-id-format       : []
 +note-id-base         : text
 +back-ref-base        : text
 +back-ref-format      : 1
 +back-ref-caret       : merge
 +notes-separator      : none
 +</refnotes>
 +~~REFNOTES~~
 +