lbteint module

This module, lbteint.py contains the integral solvers for the linearized Boltzmann Transport equations.

Contains the routines to perform the Boltzmann transport integrals.

lbteint.analytic_k_space_energy(kx, ky, kz, effmass, e_shift)

Returns the parabolic energy dispersion.

Parameters:
transport : object

A Transport() object

kx : float

The k_x in cartesian coordinates.

ky : float

The k_y in cartesian coordinates.

kz : float

The k_z in cartesian coordinates.

effmass : ndarray
Dimension: (3)

The effective mass along k_x, k_y and k_z, respectively.

Returns:
float

The energy value in eV.

Warning

This routine only accepts the diagonal elements of the effective mass tensor

lbteint.analytic_k_space_integrand(kz, ky, kx, eta, beta, effmass, e0, i, l, m)

Returns the integrand for the anlytic reciprocal space integration of the transport tensor.

Parameters:
kz : float

The k_z in cartesian coordinates.

ky : float

The k_y in cartesian coordinates.

kx : float

The k_x in cartesian coordinates.

eta : float

The reduced chemical potential.

beta : float

The \\beta factor, \\beta=(k_bT)^{-1} in eV.

effmass : float

The effective mass in units of the free electron mass.

e0 : float

The energy shift, e.g. E=\\hbar^2 k^2/2m + E_0, where E_0 is the energy shift in eV.

i : int

The order of the transport tensor.

l : {0,1,2}

The first index of the transport tensor.

m : {0,1,2}

The second index of the transport tensor.

Returns:
float

The integrand value.

lbteint.analytic_k_space_velocity(kx, ky, kz, effmass, i)

Returns the parabolic velocity dispersion.

Parameters:
kx : float

The k_x in cartesian coordinates.

ky : float

The k_y in cartesian coordinates.

kz : float

The k_z in cartesian coordinates.

effmass : ndarray
Dimension: (3)

The effective mass along k_x, k_y and k_z, respectively.

i : {0,1,2}

The direction to evaluate the velocity (0 is along k_x etc.).

Returns:
float

The velocity in eVAA.

Warning

This routine only accepts the diagonal elements of the effective mass tensor. The \\hbar/m_e factor is not returned and need to be introduced externally.

lbteint.concatenate_integrand(energies, velocities, scatter, spin_fact, chempot, beta, order)

Concatenates the integrand in the Boltzmann transport integral.

lbteint.concatenate_integrand_band(energies, velocities, tau, spin_fact, chempot, beta, order)

Concatenates the integrand in the Boltzmann transport integral and sums the bands.

lbteint.fermiintclosed(order, eta, spin_fact)

Returns the value of the closed expressions for the Fermi integrals.

lbteint.integrandpar(eps, transport, w0, eta, beta, energy_trans, effmass, i)

Returns the integrand used in the analytic energy integration of the transport coefficients in integrate_e()

Parameters:
eps : float

The reduced carrier energy.

transport : object

A Transport() object.

w0 : ndarray
Dimension: (12)

Contains the scattering rate prefactor for the different scattering mechanisms in units of inverse fs.

eta : float

The reduced chemical potential.

beta : float

The \\beta factor in eV.

energy_trans : ndarray
Dimension: (12)

Contains the energy transitions (that is added to the energy in \\tau=\\tau_0E^{r-1/2}, typically, E=E+\\hbar \\omega, where \\hbar \\omega is the size of the energy transition. Set it to zero for the non-relevant scattering mechanisms.

effmass : float

The effective mass in units of the electron mass.

i : int

The order of the transport integral to be evaluated.

Returns:
float

The integrand value.

Notes

The total scattering is calculated based on the well known scattering models for parabolic energy dispersions \\tau=\\tau_0\\epsilon^{r-1/2}, where r is the scattering factor.

lbteint.integrandpardos(eps, transport, w0, eta, beta, energy_trans, effmass, i)

The integrand for the density of states integral over energy.

Parameters:
eps : float

The reduced carrier energy

transport : object

A Transport() object

w0 : ndarray
Dimension: (12)

Contains the scattering rate prefactor for the different scattering mechanisms in units of inverse fs. Not used in this routine, but it needs the dummy from the call argument.

eta : float

The reduced chemical potential

beta : float

The \\beta factor in eV. Not used in this routine, but it needs the dummy from the call argument.

energy_trans : ndarray
Dimension: (12)

Contains the energy transitions (that is added to the energy in \\tau=\\tau_0E^{r-1/2}, typically, E=E+\\hbar \\omega, where \\hbar \\omega is the size of the energy transition. Set it to zero for the non-relevant scattering mechanisms. Not used in this routine, but it needs the dummy from the call argument.

effmass : float

The effective mass in units of the electron mass. Not used in this routine, but it needs the dummy from the call argument.

i : int

The order of the transport integral to be evaluated. Not used in this routine, but it needs the dummy from the call argument.

Returns:
float

The integrand value for the density of states.

Notes

Calculates the density of states integrand

System Message: WARNING/2 (\\frac{\\epsilon^{1/2}}{1+\\exp{(\\epsilon-\\eta)}} )

latex exited with error [stdout] This is pdfTeX, Version 3.14159265-2.6-1.40.18 (TeX Live 2017/Debian) (preloaded format=latex) restricted \write18 enabled. entering extended mode (./math.tex LaTeX2e <2017-04-15> Babel <3.18> and hyphenation patterns for 84 language(s) loaded. (/usr/share/texlive/texmf-dist/tex/latex/base/article.cls Document Class: article 2014/09/29 v1.4h Standard LaTeX document class (/usr/share/texlive/texmf-dist/tex/latex/base/size12.clo)) (/usr/share/texlive/texmf-dist/tex/latex/base/inputenc.sty (/usr/share/texlive/texmf-dist/tex/latex/ucs/utf8x.def)) (/usr/share/texlive/texmf-dist/tex/latex/ucs/ucs.sty (/usr/share/texlive/texmf-dist/tex/latex/ucs/data/uni-global.def)) (/usr/share/texlive/texmf-dist/tex/latex/amsmath/amsmath.sty For additional information on amsmath, use the `?' option. (/usr/share/texlive/texmf-dist/tex/latex/amsmath/amstext.sty (/usr/share/texlive/texmf-dist/tex/latex/amsmath/amsgen.sty)) (/usr/share/texlive/texmf-dist/tex/latex/amsmath/amsbsy.sty) (/usr/share/texlive/texmf-dist/tex/latex/amsmath/amsopn.sty)) (/usr/share/texlive/texmf-dist/tex/latex/amscls/amsthm.sty) (/usr/share/texlive/texmf-dist/tex/latex/amsfonts/amssymb.sty (/usr/share/texlive/texmf-dist/tex/latex/amsfonts/amsfonts.sty)) (/usr/share/texlive/texmf-dist/tex/latex/anyfontsize/anyfontsize.sty) (/usr/share/texlive/texmf-dist/tex/latex/tools/bm.sty) (./math.aux) (/usr/share/texlive/texmf-dist/tex/latex/ucs/ucsencs.def) (/usr/share/texlive/texmf-dist/tex/latex/amsfonts/umsa.fd) (/usr/share/texlive/texmf-dist/tex/latex/amsfonts/umsb.fd) ! Missing } inserted. <inserted text> } l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Missing } inserted. <inserted text> } l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Extra }, or forgotten $. <argument> ...n^{1/2}}{1+\\exp{(\\epsilon-\\eta)}} \end {split} l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Missing { inserted. <inserted text> { l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Missing } inserted. <inserted text> } l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Missing } inserted. <inserted text> } l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Extra }, or forgotten $. <argument> ...n^{1/2}}{1+\\exp{(\\epsilon-\\eta)}} \end {split} l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} ! Missing { inserted. <inserted text> { l.14 ...2}}{1+\\exp{(\\epsilon-\\eta)}}\end{split} [1] (./math.aux) ) (see the transcript file for additional information) Output written on math.dvi (1 page, 532 bytes). Transcript written on math.log.
lbteint.integrandpart2(eps, transport, w0, eta, beta, energy_trans, effmass, i)

Returns the integrand used in the analytic energy integration of the transport distribution function with a quadratic \\tau term

Parameters:
eps : float

The reduced carrier energy.

transport : object

A Transport() object.

w0 : ndarray
Dimension: (12)

Contains the scattering rate prefactor for the different scattering mechanisms in units of inverse fs.

eta : float

The reduced chemical potential.

beta : float

The \\beta factor in eV.

energy_trans : ndarray
Dimension: (12)

Contains the energy transitions (that is added to the energy in \\tau=\\tau_0E^{r-1/2}, typically, E=E+\\hbar \\omega, where \\hbar \\omega is the size of the energy transition. Set it to zero for the non-relevant scattering mechanisms.

effmass : float

The effective mass in units of the electron mass.

i : int

The order of the transport integral to be evaluated.

Returns:
float

The integrand value.

Notes

The total scattering is calculated based on the well known scattering models for parabolic energy dispersions \\tau=\\tau_0\\epsilon^{r-1/2}, where r is the scattering factor. Difference from integrandpar(): here tau^2 is used in the integrand (for the calculation of the Hall factor).

lbteint.scipy_e_integrals(transport, integrand, e_min, e_max, w0, eta, beta, energy_trans, effmass, order, spin_fact, method='quad')

Calculates the one dimensional energy integrals.

Uses the SciPy function scipy.integrate.quad().

Parameters:
transport : object

A Transport() object

integrand : {“normal”,”hall”,”dos”}

Selects the type of integrand to be used. “normal” selects integrandpar(). “hall” selects integrandpart2(), “dos” selects integrandpardos().

e_min : float

The lower integration limit in eV.

e_max : float

The higher integration limit in eV.

w0 : ndarray
Dimension: (12)

Contains the scattering rate prefactor (inverse of relaxation time) for the different scattering mechanisms in units of inverse fs.

eta : float

The reduced chemical potential.

beta : float

The \\beta factor, (\\mathrm{k_b}T)^{-1} in eV.

effmass : ndarray
Dimension: (3)

The effective mass along the three reciprocal unit vectors in units of the free electron mass.

e0 : float

The energy shift, e.g. E=\\hbar^2 k^2/2m + E_0, where E_0 is the energy shift in eV.

i : int

The order of the transport tensor.

l : {0,1,2}

The first index of the transport tensor.

m : {0,1,2}

The second index of the transport tensor.

spin_fact : int

The spin degeneracy factor. 1 for non-spin degeneracy, 2 for spin degeneracy.

method : {“quad”}, optional

The SciPy three dimensional integration method using scipy.integrate.quad().

Returns:
float

The resulting integral over energy.

lbteint.scipy_k_integrals(eta, beta, effmass, e0, i, l, m, method='tplquad')

Calculates the three dimensional wave vector integrals.

Uses the SciPy function scipy.integrate.tplquad().

Parameters:
eta : float

The reduced chemical potential

beta : float

The \\beta factor, (\\mathrm{k_b}T)^{-1} in eV.

effmass : ndarray
Dimension: (3)

The effective mass along the three reciprocal unit vectors in units of the free electron mass.

e0 : float

The energy shift, e.g. E=\\hbar^2 k^2/2m + E_0, where E_0 is the energy shift in eV.

i : int

The order of the transport tensor.

l : {0,1,2}

The first index of the transport tensor

m : {0,1,2}

The second index of the transport tensor

method : {“tplquad”}, optional

The SciPy three dimensional integration method using scipy.integrate.tplquad().

Returns:
float

The resulting integral over the wave vectors.

lbteint.scipy_k_integrals_discrete(tr, integrand_type, energies, velocities, scatter, chempot, beta, order, spin_fact, method='trapz')

Calculates the three dimensional integrals over the k-points for discrete data.

Uses SciPy integration functions for discrete data.

Parameters:
tr : object

A Transport() object

energies: ndarray

Contains the band energies in eV for each k-point.

velocities: ndarray

Contains the derivative if energies without the \hbar factors for each k-point.

scatter:

Contains the relaxation time at each k-point.

chempot : float

The chemical potential in eV

beta : float

The \\beta factor, (\\mathrm{k_b}T)^{-1} in eV.

spin_fact : int

The spin factor, 1 for non-spin degeneracy and 2 for spin degeneracy.

method : {“trapz”, “simps”, “romb”}, optional

The SciPy three dimensional integration method for the scipy.integrate.trapz(), scipy.integrate.simps() and the scipy.integrate.romb() functions, respectively. Defaults to “trapz”.

Returns:
integral : float

The resulting integral over the wave vectors.

lbteint.scipy_k_integrals_discrete2(tr, energies, velocities, scatter, chempot, beta, spin_fact, order, method='trapz')

Calculates the three dimensional integrals over the k-points for discrete data.

Uses integration functions for discrete data.

Parameters:
tr : object

A Transport() object

chempot : float

The chemical potential in eV

beta : float

The \\beta factor, (\\mathrm{k_b}T)^{-1} in eV.

spin_fact : int

The spin factor, 1 for non-spin degeneracy and 2 for spin degeneracy.

kx, ky, kz : float, float, float

The spacing in inverse AA between the points along each direction.

order : float

The order of the energy minus chemical potential term in the denominator.

method : {“trapz”, “simps”, “romb”}, optional

The SciPy three dimensional integration method for the scipy.integrate.trapz(), scipy.integrate.simps() and the scipy.integrate.romb() functions, respectively. Defaults to “trapz”.