atoMEC.postprocess package

Submodules

atoMEC.postprocess.conductivity module

The conductivity module handles routines used to model the electrical conductivity.

So far just the Kubo-Greenwood method is implemented.

Classes

  • KuboGreenwood : Holds various routines needed to compute the Kubo-Greenwood conducivity, including its various components. Also contains various properties related to the KG conducivity.

  • SphHamInts: Holds the routines to construct angular integrals from the spherical harmonic functions

  • RadialInts: Holds the routines to construct radial integrals from the radial KS orbitals

class atoMEC.postprocess.conductivity.KuboGreenwood(Atom, model, orbitals, valence_orbs=[], nmax=0, lmax=0)[source]

Bases: object

Class for Kubo-Greenwood conductivity and MIS via TRK sum rule.

static calc_eig_diff_mat(eigvals, orb_subset_1, orb_subset_2)[source]

Compute the matrix of eigenvalue differences e_l1n1 - e_ln2n2.

Parameters:
  • eigvals (ndarray) – the KS energy eigenvalues

  • orb_subset_1 (tuple) – the first subset of orbitals (eg valence)

  • orb_subset_2 (tuple) – the second subset of orbitals (eg conduction)

Returns:

occ_diff_mat – the occupation number difference matrix

Return type:

ndarray

static calc_mel_grad_int(orb_l1n1, orb_l2n2, l1, n1, l2, n2, m, xgrid)[source]

Calculate the matrix element \(|<\phi_{n1l1}|\nabla|\phi_{n1l2}>|^2\).

Parameters:
  • orb_l1n1 (ndarray) – l1,n1 radial KS orbital

  • orb_l2n2 (ndarray) – l2,n2 radial KS orbital

  • l1 (int) – 1st angular momentum quantum number

  • n1 (int) – 1st principal quantum number

  • l2 (int) – 2nd angular momentum quantum number

  • n2 (int) – 2nd principal quantum number

  • m (int) – magnetic quantum number

  • xgrid (ndarray) – log grid

Returns:

mel_grad_int – the matrix element \(|<\phi_{n1l1}|\nabla|\phi_{n1l2}>|^2\).

Return type:

float

static calc_occ_diff_mat(occnums, orb_subset_1, orb_subset_2)[source]

Compute the matrix of occupation number diffs -(f_l1n1 - f_l2n2).

Parameters:
  • occnums (ndarray) – the (unweighted FD) KS occupation numbers

  • orb_subset_1 (tuple) – the first subset of orbitals (eg valence)

  • orb_subset_2 (tuple) – the second subset of orbitals (eg conduction)

Returns:

occ_diff_mat – the occupation number difference matrix

Return type:

ndarray

calc_sig(R1_int, R2_int, orb_subset_1, orb_subset_2)[source]

Compute the integrated dynamical conducivity for given subsets (see notes).

Parameters:
  • R1_int (ndarray) – the ‘R1’ radial component of the integrand (see notes)

  • R2_int (ndarray) – the ‘R2’ radial component of the integrand (see notes)

  • orb_subset_1 (list of tuples) – the first subset of orbitals to sum over

  • orb_subset_2 (list of tuples) – the second subset of orbitals to sum over

Returns:

sig – the integrated dynamical conductivity

Return type:

float

Notes

This function returns the integrated dynamical conductivity, \(\bar{\sigma}=\int_0^\infty d\omega \sigma(\omega)\). The conductivity \(\sigma(\omega)\) is defined as

\[\sigma_{S_1,S2}(\omega) = \frac{2\pi}{3V\omega} \sum_{i\in S_1}\sum_{j\in S_2} (f_i - f_j)\ |\langle\phi_{i}|\nabla|\phi_{j}\rangle|^2\delta(\epsilon_j-\epsilon_i-\omega),\]

where \(S_1,S_2\) denote the subsets of orbitals specified in the function’s paramaters, e.g. the conduction-conduction orbitals.

In practise, the integral in the above equation is given by a discrete sum due to the presenence of the dirac-delta function.

The paramaters R1_int and R2_int refer to radial integral components in the calculation of the matrix elements. See the supplementary information of Ref. [8] for more information on thse components, and the functions calc_R1_int_mat() and calc_R2_int_mat() for their definitions.

References

calc_sig_func(R1_int, R2_int, orb_subset_1, orb_subset_2, omega_max, n_freq, gamma)[source]

Compute the dynamical conducivity for given subsets (see notes).

Parameters:
  • R1_int (ndarray) – the ‘R1’ radial component of the integrand (see notes)

  • R2_int (ndarray) – the ‘R2’ radial component of the integrand (see notes)

  • orb_subset_1 (list of tuples) – the first subset of orbitals to sum over

  • orb_subset_2 (list of tuples) – the second subset of orbitals to sum over

  • omega_max (float) – maximum value of the frequency grid

  • n_freq (int) – number of points in the frequency grid

  • gamma (float) – smoothing factor for the Lorentzian

Returns:

sig_omega, nele

  • sig_omega: 2d array containing frequency grid and conductivity \(\sigma(\omega)\)

  • n_ele: the number of electrons from integration of \(\sigma(\omega)\); equivalent to N_ij (for orb subsets ij) in the limit \(\gamma\to 0\)

Return type:

tuple (ndarray, float)

Notes

This function returns the dynamical conductivity, \(\sigma(\omega)\), defined as

\[\begin{split}\sigma_{S_1,S2}(\omega) &= \frac{2\pi}{3V\omega} \sum_{i\in S_1}\sum_{j\in S_2} (f_i - f_j)\ |\langle\phi_{i}|\nabla|\phi_{j}\rangle|^2\ \mathcal{L}(\epsilon_i, \epsilon_j, \gamma, \omega) \\ \mathcal{L}(\epsilon_i, \epsilon_j, \gamma, \omega)&=\ \frac{\gamma}{\pi}\frac{1}{\gamma^2+(\omega+[\epsilon_i-\epsilon_j)])^2}\end{split}\]

where \(S_1,S_2\) denote the subsets of orbitals specified in the function’s paramaters, e.g. the conduction-conduction orbitals.

As can be seen in the above equation, the dirac-delta function in the definition of the KG conductivity (see calc_sig function) is represented by a Lorentzian distribution \(\mathcal{L}\) to obtain a smooth conductivity function. In the limit \(\gamma\to 0\), the Lorentzian becomes a delta function.

The paramaters R1_int and R2_int refer to radial integral components in the calculation of the matrix elements. See the supplementary information of Ref. [8] for more information on these components, and the functions calc_R1_int_mat() and calc_R2_int_mat() for their definitions.

check_sum_rule(l, n, m)[source]

Check the sum rule (see notes) for an orbital \(\phi_{nlm}\) is satisfied.

Parameters:
  • l (int) – angular quantum number

  • n (int) – principal quantum number

  • m (int) – magnetic quantum number

Returns:

sum_mom – the momentum sum rule (see notes)

Return type:

ndarray

Notes

The expression for the momentum sum rule is given by

\[S_{p} = \sum_{(n_1,l_1,m_1)\neq (n,l,m)}\ \frac{|\langle\phi_{nlm}|\nabla|\phi_{n_1 l_1 m_1}\rangle|^2} {\ \epsilon_{n_1,l_1,m_1}-\epsilon_{n,l,m}}\]

If the sum rule is satisfied, the summation above should equal 1/2. See Eq. (38) of Ref. [7] for an explanation of this sum rule.

References

cond_tot(component='tt', gamma=0.01, maxfreq=50, nfreq=500)[source]

Calculate the chosen component of dynamical electrical conductivity sig(w).

Parameters:
  • component (str, optional) – the desired component of the conducivity e.g. “cc”, “tt” etc

  • gamma (float, optional) – smoothing factor

  • maxfreq (float, optional) – maximum frequency to scan up to

  • nfreq (int, optional) – number of points in the frequency grid

Returns:

cond_tot_ – dynamical electrical conductivity

Return type:

ndarray

static sig_to_N(sig, V)[source]

Map the integrated conducivity to electron number.

Parameters:
  • sig (float) – integrated conducivity

  • V (float) – volume of sphere

Returns:

N_ele – electron number

Return type:

float

property N_free

the free electron number from TRK sum-rule.

Type:

float

property N_tot

the total electron number from TRK sum-rule.

Type:

float

property R1_int_cc

Conducting-conducting component of the R1 radial integral.

property R1_int_cv

Conducting-valence component of the R1 radial integral.

property R1_int_tt

Total-total component of the R1 radial integral.

property R1_int_vv

Valence-valence component of the R1 radial integral.

property R2_int_cc

Conducting-conducting component of the R2 radial integral.

property R2_int_cv

Conducting-valence component of the R2 radial integral.

property R2_int_tt

Total-total component of the R2 radial integral.

property R2_int_vv

Valence-valence component of the R2 radial integral.

property all_orbs

all the possible orbital pairings.

Type:

List of tuples

property cond_orbs

all the conduction band orbital pairings.

Type:

List of tuples

property sig_cc

the integrated cc conductivity component.

Type:

ndarray

property sig_cv

the integrated cv conductivity component.

Type:

ndarray

property sig_tot

the integrated total conductivity.

Type:

ndarray

property sig_vv

the integrated vv conductivity component.

Type:

ndarray

property sph_vol

the volume of the sphere.

Type:

float

class atoMEC.postprocess.conductivity.RadialInts[source]

Bases: object

Contains functions required to compute various integrals of the radial KS fns.

static R1_int_term(eigfunc, grad_orb2, xgrid)[source]

Input function to the calc_R1_int_mat() function.

Parameters:
  • eigfunc (ndarray) – KS orbital l1,n1

  • grad_orb2 (ndarray) – derivative of KS orbital l2,n2

  • xgrid (ndarray) – log grid

Returns:

R1_int – the matrix element for the R1_int_mat function

Return type:

float

static R2_int_term(eigfunc_1, eigfunc_2, xgrid)[source]

Input function to the calc_R2_int_mat() function.

Parameters:
  • eigfunc_1 (ndarray) – KS orbital l1,n1

  • eigfunc_2 (ndarray) – KS orbital l2,n2

  • xgrid (ndarray) – log grid

Returns:

R2_int – the matrix element for the R2_int_mat function

Return type:

float

static calc_R1_int(orb1, orb2, xgrid)[source]

Compute the R1 integral between two orbitals orb1 and orb2 (see notes).

Parameters:
  • orb1 (ndarray) – the first radial orbital

  • orb2 (ndarray) – the second radial orbital

Returns:

R1_int – the R1 integral

Return type:

ndarray

Notes

See calc_R1_int_mat() for definition of the integral

classmethod calc_R1_int_mat(eigfuncs, occnums, xgrid, orb_subset_1, orb_subset_2)[source]

Compute the ‘R1’ integral matrix (see notes).

Parameters:
  • eigfuncs (ndarray) – the KS eigenfunctions

  • occnums (ndarray) – the KS occupation numbers

  • xgrid (ndarray) – the log grid

  • orb_subset_1 (tuple) – the first subset of orbitals (eg valence)

  • orb_subset_2 (tuple) – the second subset of orbitals (eg conduction)

Returns:

R1_mat – the R1 integral matrix (see notes)

Return type:

ndarray

Notes

The definition of the R1 integral is (see Ref. [7] and supplementary of [8])

\[R^{(1)}=4\pi\int_0^R dr r^2 X_{n_1 l_1}(r) \frac{dX_{n_2 l_2}(r)}{dr},\]

where \(X_{nl}(r)\) are the radial KS functions.

static calc_R2_int(orb1, orb2, xgrid)[source]

Compute the R2 integral between two orbitals orb1 and orb2 (see notes).

Parameters:
  • orb1 (ndarray) – the first radial orbital

  • orb2 (ndarray) – the second radial orbital

Returns:

R2_int – the R2 integral

Return type:

ndarray

Notes

See calc_R2_int_mat() for definition of the integral

classmethod calc_R2_int_mat(eigfuncs, occnums, xgrid, orb_subset_1, orb_subset_2)[source]

Compute the ‘R2’ integral matrix (see notes).

Parameters:
  • eigfuncs (ndarray) – the KS eigenfunctions

  • occnums (ndarray) – the KS occupation numbers

  • xgrid (ndarray) – the log grid

  • orb_subset_1 (tuple) – the first subset of orbitals (eg valence)

  • orb_subset_2 (tuple) – the second subset of orbitals (eg conduction)

Returns:

R2_mat – the R2 integral matrix (see notes)

Return type:

ndarray

Notes

The definition of the R2 integral is (see Ref. [7] and supplementary of [8])

\[R^{(1)}=4\pi\int_0^R dr r X_{n_1 l_1}(r) X_{n_2 l_2}(r),\]

where \(X_{nl}(r)\) are the radial KS functions.

class atoMEC.postprocess.conductivity.SphHamInts[source]

Bases: object

Contains the functions needed to compute various spherical harmonic integrals.

static P2_func(x, l1, l2, m)[source]

Calculate the ‘P2’ function (see notes).

Parameters:
  • x (float) – input for Legendre polynomial

  • l1 (int) – 1st angular quantum number

  • l2 (int) – 2nd angular quantum number

  • m (int) – magnetic quantum number

Returns:

P2_func_ – the P2 function

Return type:

float

Notes

The P2 function is defined as (see also Refs. [7] and [8])

\[f_p^{(2)}[l_1,l_2,m](x) = x P_{l_1}^m (x) P_{l_2}^m (x)\]

where P_{l}^m (x) are Legendre polynomial functions.

static P4_func(x, l1, l2, m)[source]

Calculate the ‘P4’ function (see notes).

Parameters:
  • x (float) – input for Legendre polynomial

  • l1 (int) – 1st angular quantum number

  • l2 (int) – 2nd angular quantum number

  • m (int) – magnetic quantum number

Returns:

P4_func_ – the P4 function

Return type:

float

Notes

The P4 function is defined as (see also Refs. [7] and [8])

\[\begin{split}f_p^{(4)}[l_1,l_2,m](x)&=-(1-x)^2 P^m_{l_1}(x) \frac{dP_{l_2}^m(x)}{dx}\\ &= P^m_{l_1}(x) [(l_2+m)P_{l_2-1}^m(x)-xl_2\ P_{l_2}^m(x)]\end{split}\]

where \(P_{l}^m(x)\) are Legendre polynomial functions.

classmethod P_int(func_int, l1, l2, m)[source]

Integrate the P2 or P4 function (see notes).

Parameters:
  • func_int (int) – the desired P integral (can be 2 or 4)

  • l1 (int) – 1st angular quantum number

  • l2 (int) – 2nd angular quantum number

  • m (int) – magnetic quantum number

Returns:

P_int_ – the integrated P2 or P4 function

Return type:

float

Notes

The integrals are defined as

\[\bar{P}^{(n)}_{ll'm} = 2\pi c_{lm}c_{l'm}\int_{-1}^1 dx \ f_p^{(n)}[l_1,l_2,m](x)\]

With the functions \(f_p^{(n)}(x)\) defined below (P2_func() and P4_func()).

classmethod P_mat_int(func_int, lmax)[source]

Compute the matrix of P function (angular) integrals (see notes).

Parameters:
  • func_int (int) – the desired P integral (can be 2 or 4)

  • lmax (int) – the maximum value of angular momentum

Returns:

P_mat – matrix of P func integrals for chosen func_int

Return type:

ndarray

Notes

See Refs. [7] and [8] (supplemental material) for the definitions of the P2 and P4 functions, ands the P2_func(), P4_func() and P_int() functions.

static sph_ham_coeff(l, m)[source]

Compute coefficients of spherical harmonic functions.

Parameters:
  • l (int) – angular quantum number

  • m (int) – magnetic quantum number

Returns:

c_lm – coefficient for spherical harmonic function (l,m) (see notes)

Return type:

float

Notes

The spherical harmonic functions with coefficients \(c_{lm}\) are defined as

\[\begin{split}Y_m^l(\theta,\phi) &= c_{lm} P_l^m (\cos\theta) e^{im\phi}\\ c_{lm} &= \sqrt{\frac{(2l+1)(l-m)!}{4\pi(l+m)!}}\end{split}\]

atoMEC.postprocess.localization module

The localization module handles routines designed to measure electron localization.

Classes

  • ELFTools : Holds functions to compute the electron localization function (ELF) and related quantities e.g. number of electrons per shell.

Functions

  • calc_IPR_mat(): Computes the inverse participation ratio (IPR) for the orbitals.

class atoMEC.postprocess.localization.ELFTools(Atom, model, orbitals, density, method='orbitals')[source]

Bases: object

ELF is the electron localization function, a measure of electron localization.

It can be used as a tool for qualitative insight as well as an additional method to compute the mean ionization state. This class contains routines to compute the ELF and other related useful properties.

Parameters:
static calc_ELF(eigfuncs, occnums, xgrid, density, method='orbitals')[source]

Compute the ELF (see notes).

Parameters:
  • eigfuncs (ndarray) – the (modified) KS eigenfunctions \(P(x)=\exp(x/2)X(x)\)

  • occnums (ndarray) – the orbital occupations

  • xgrid (ndarray) – the logarithmic grid

  • density (ndarray) – the electron density

Returns:

ELF – the electron localization function

Return type:

ndarray

Notes

The ELF is defined as

\[\begin{split}\mathrm{ELF}^\sigma(r) = \frac{1}{1 + (\chi^\sigma(r))^2}\,, \\ \chi^\sigma(r) = D^\sigma(r) / D^\sigma_0(r),\end{split}\]

where \(D^\sigma(r)\) and \(D^\sigma_0(r)\) are the electron pair density functions for the system of interest and the uniform electron gas (UEG) respectively.

static calc_N_shell(xargs_min, xgrid, density, ELF)[source]

Compute the number of electrons in each shell from the ELF.

The number of electrons per shell is determined by integrating the electron density between successive minima of the ELF.

Parameters:
  • xargs_min (list) – list of minima, as well as the first and last points

  • xgrid (ndarray) – the logarithmic grid

  • density (ndarray) – the electron density

  • ELF (ndarray) – the electron localization function

Returns:

N_shell – list of the number of electrons per shell

Return type:

list of floats

static calc_epdc(eigfuncs, occnums, xgrid, density)[source]

Calculate the electron pair density curvature (see notes).

Parameters:
  • eigfuncs (ndarray) – the (modified) KS eigenfunctions \(P(x)=\exp(x/2)X(x)\)

  • occnums (ndarray) – the orbital occupations

  • xgrid (ndarray) – the logarithmic grid

  • density (ndarray) – the electron density

Returns:

epdc – the electron pair density curvature

Return type:

ndarray

Notes

The epdc is defined as

\[D^\sigma(r) = \tau^\sigma(r) - \frac{1}{8}\frac{|\nabla\rho(r)|^2}{\rho(r)},\]

where \(\tau^\sigma(r)\) is the local kinetic energy density.

static calc_epdc_dens(xgrid, density)[source]

Calculate the electron pair density curvature (see notes).

Parameters:
  • xgrid (ndarray) – the logarithmic grid

  • density (ndarray) – the electron density

Returns:

epdc – the electron pair density curvature

Return type:

ndarray

Notes

The epdc is defined as

\[D^\sigma(r) = D_0^\sigma(r) - \frac{1}{9}\frac{|\nabla\rho(r)|^2}{\rho(r)}\ + \frac{1}{6} \nabla^2 \rho(r),\]

where \(D_0^\sigma(r)=(3/10)(3\pi^2)^{2/3} \rho^{5/3}(r)\) is the epdc for the UEG.

static get_minima(ELF, xgrid)[source]

Get the locations of the minima of the ELF.

Parameters:
  • ELF (ndarray) – the electron localization function

  • xgrid (ndarray) – the logarithmic grid

Returns:

xargs_min – list of locations of the minima and first and last points

Return type:

list of ints

property ELF

the electron localization function for a given system.

Type:

ndarray

property N_shell

the number of electrons in each ‘shell’.

Type:

ndarray

property epdc

the electron pair density curvature.

Type:

ndarray

atoMEC.postprocess.localization.MIS_count(model, orbitals, core_orbs)[source]

Calculate the MIS using the “counting” method.

Parameters:
  • model (models.ISModel) – the ISModel object

  • orbitals (staticKS.orbitals) – the KS orbitals object

  • core_orbs (list of tuples) – the core orbitals

Returns:

MIS – the mean ionization state

Return type:

np.ndarray

atoMEC.postprocess.localization.calc_IPR_mat(eigfuncs, xgrid, grid_type=None)[source]

Calculate the inverse participation ratio for all eigenfunctions (see notes).

Parameters:
  • eigfuncs (ndarray) – transformed radial KS orbitals \(P_{nl}(x)=\exp(x/2)X_{nl}(x)\)

  • xgrid (ndarray) – the logarithmic grid

Returns:

IPR_mat – the matrix of all IPR values

Return type:

ndarray

Notes

The inverse participation ratio for a state i is usually defined as

\[\mathrm{IPR}_i = \int \mathrm{d}{\mathbf{r}} |\Psi_i(\mathbf{r})|^4\]

It is typically used as a localization measure.

Warning

The current definition in this version is not mathematically correct. It does not include the proper contribution from the spherical harmonics \(|Y_l^m(\theta,\phi)|^4\). This is omitted as it makes little difference to the flavour of the results but complicates things. Currently, this function does not correctly produce the expected limits (even if the spherical harmonic contribution is correctly accounted for). Use at your own peril…

atoMEC.postprocess.pressure module

The pressure module contains functions to compute the pressure with various approaches.

Functions

  • finite_diff() : Calculate electronic pressure with finite-difference method.

  • stress_tensor() : Calculate electronic pressure with stress-tensor method.

  • virial() : Calculate electronic pressure with virial method.

  • ideal_electron() : Calculate electronic pressure with ideal method.

  • calc_Wd_xc() : Calculate the derivative contribution to virial pressure.

  • ions_ideal() : Calculate ionic pressure with ideal gas law.

atoMEC.postprocess.pressure.calc_Wd_xc(xc_func_id, density, grid_type)[source]

Compute the ‘derivative’ component of the xc term in virial formula (see notes).

Parameters:
  • xc_func_id (str) – the exchange or correlation functional id

  • density (staticKS.Density) – the density object

Returns:

Wd_xc – the derivative component of the xc term in the virial formula.

Return type:

float

Notes

In the LDA, Wd_xc is given by [10]

\[W^\mathrm{d}_\mathrm{xc} = \int \mathrm{d}r v_\mathrm{xc}(r) n(r)\]

GGA implementation is still to come.

atoMEC.postprocess.pressure.finite_diff(atom, model, orbs, pot, conv_params={}, scf_params={}, force_bound=[], write_info=False, verbosity=0, dR=0.01, method='A')[source]

Calculate the electronic pressure using the finite differences method.

Parameters:
  • atom (atoMEC.Atom) – The main atom object

  • model (models.ISModel) – The ISModel object

  • orbs (staticKS.Orbitals) – the orbitals object

  • pot (staticKS.Potential) – the potential object

  • conv_params (dict, optional) – dictionary of convergence parameters as follows: { econv (float) : convergence for total energy, nconv (float) : convergence for density, vconv (float) : convergence for electron number, eigtol (float) : tolerance for eigenvalues }

  • scf_params (dict, optional) – dictionary for scf cycle parameters as follows: { maxscf (int) : maximum number of scf cycles, mixfrac (float) : density mixing fraction }

  • force_bound (list of list of ints, optional) – force certain levels to be bound, for example: force_bound = [0, 1, 0] forces the orbital quith quantum numbers \(\sigma=0,\ l=1,\ n=0\) to be always bound even if it has positive energy. This prevents convergence issues.

  • verbosity (int, optional) – how much information is printed at each SCF cycle. verbosity=0 prints the total energy and convergence values (default) verbosity=1 prints the above and the KS eigenvalues and occupations.

  • write_info (bool, optional) – prints the scf cycle and final parameters

  • dR (float, optional) – radius difference for finite difference calculation

  • method (str, optional) – method for computing the free energy: can either use normal construction (“A”) or with the EnergyAlt class (“B”)

Returns:

P_e – electronic pressure in Ha

Return type:

float

atoMEC.postprocess.pressure.ideal_electron(Atom, chem_pot)[source]

Compute the ideal electron pressure.

Parameters:
  • Atom (atoMEC.Atom) – the atom object

  • chem_pot (float) – the chemical potential

Returns:

P_e – the ideal electron pressure

Return type:

float

Notes

The formula to determine the ideal electron pressure is

\[P_\textrm{e} = \frac{2^{3/2}}{3\pi^2} \int\mathrm{d}\epsilon \epsilon^{3/2} f_\textrm{FD}(\epsilon,\beta,\mu)\]
atoMEC.postprocess.pressure.ions_ideal(atom)[source]

Compute the the ideal gas pressure for the ions.

In atomic units, the ideal gas pressure is simply P = T/V.

Parameters:

atom (atoMEC.Atom) – the Atom object

Returns:

P_ion – the ionic pressure

Return type:

float

atoMEC.postprocess.pressure.stress_tensor(Atom, model, orbs, pot, only_rr=False)[source]

Calculate the pressure with the stress tensor approach [9].

Parameters:
  • Atom (atoMEC.Atom) – The main atom object

  • model (models.ISModel) – The model object

  • orbs (staticKS.Orbitals) – the orbitals object

  • pot (staticKS.Potential) – the potential object

  • only_rr (bool, optional) – whether to use just the radial component of the stress tensor (True) or the full trace (False). See [9] for definitions.

Returns:

P_e – the electronic pressure

Return type:

float

References

atoMEC.postprocess.pressure.virial(atom, model, energy, density, orbs, pot, use_correction=False, method='A')[source]

Compute the pressure using the virial theorem (see notes).

Parameters:
Returns:

P_e – the electronic pressure

Return type:

float

Notes

The virial pressure is given by the formula [10] [12]

\[\begin{split}P &= \frac{K1 + K2 + E_\mathrm{en} + E_\mathrm{Ha} + W_\mathrm{xc}}{3V}\ , \\ W_\mathrm{xc} &= 3 (W^\mathrm{d}_\mathrm{xc} - E_\mathrm{xc})\end{split}\]

If use_correction==True, \(K1\) and \(K2\) are respectively the integrated kinetic energy densities “B” and “A” as in staticKS.Energy.calc_E_kin_dens().

If use_correction==False, both terms \(K1\) and \(K2\) are the same and both given by method “A” in staticKS.Energy.calc_E_kin_dens().

References

Module contents

Everything that depends on the output from an SCF calculation.