Documentation
Types
Mstar2t.Utils.ParabBand
— TypeParabBand(mstar::Vector{Float64},ϵ₀::Union{Float64,Vector{Float64}},type::Int64,deg::Int64)
Julia type to represent a band. Parameters:
- mstar: effective mass tensor
- ϵ₀: minimum/maximum of the band
- type: cconduction +1 or valence -1 band
- deg: degeneracy
Example
julia> mstar = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0];
julia> ϵ₀ = 0.0;
julia> type = 1;
julia> deg = 1;
julia> band = ParabBand(mstar,ϵ₀,type,deg); # band creation
Mstar2t.Utils.BandStructure
— TypeBandStructure(n::Int64, bands::Vector{ParabBand}, μ::Union{Float64,Vector{Float64}})
Julia type to represent a band full band structure. Parameters:
- n: number of bands
- bands: vector of bands
- μ: Fermi level position
Example
julia> m_1 = [1., 2., 3., 0.0, 0.0, 0.0];
julia> m_2 = [.5, .7, .9, 0.0, 0.0, 0.0];
julia> band_1 = ParabBand(m_1,0.0,1,1); # create the conduction band
julia> band_2 = ParabBand(m_2,-0.1,-1,1);
julia> μ = collect(-.1:.01:.1);
julia> model = BandStructure(2,[band_1,band_2],μ); # band creation
Functions
Mstar2t.Transport.electrical_conductivity
— Functionelectrical_conductivity(bandstructure::BandStructure, Ts::Union{Vector{Float64},Float64,Vector{Int64},Int64}, τ::Union{ScModel,Matthiessen}; exportasdf::Bool=false, fulltensor::Bool=false)
Function that by default computes the trace of the tensorial version of the electical conductivity in units of $(\Omega m)^{-1}$ for a given bandstructure
and a choice of the relaxation time $\tau$, for given values of chemical potential $\mu$ and temperature Ts
. The transport coefficient is returned as a matrix
of dimensions (length(T),length(μ))
. The boolean variable exportasdf
allows to return the calculations as a DataFrame
with all the parameters included. If fulltensor
is set to true
the full tensor is returned in place of the trace. Reference: "Theory of band warping and its effects on thermoelectronic transport properties", PHYSICAL REVIEW B89.
Example
julia> m = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0]; # effectuve mass
julia> band = ParabBand(m,0.0,1,1); # build the band
julia> μ = collect(-.1:.01:.1); # fermi level position
julia> model = BandStructure(1,band_1,μ) # build the band structure
julia> T = collect(50:10:650); # temperature range
julia> τ_form = Scattering.constant() # relaxation time
julia> σ = electrical_conductivity(model,T,τ_form); # compute the electrical conductivity
Mstar2t.Transport.seebeck_coefficient
— Functionseebeck_coefficient(bandstructure::BandStructure, Ts::Union{Vector{Float64},Float64,Vector{Int64},Int64},τ::Union{ScModel,Matthiessen}; exportasdf::Bool=false, fulltensor::Bool=false)
Function that by default computes the Seebeck coefficient in units of $V/K$ for a given bandstructure
and a choice of the relaxation time $\tau$, for given values of chemical potential $\mu$ and temperature Ts
. The transport coefficient is returned as a matrix
of dimensions (length(T),length(μ))
. The boolean variable exportasdf
allows to return the calculations as a DataFrame
with all the parameters included. If fulltensor
is set to true
, the full tensor is returned in place of the trace. Reference: "Theory of band warping and its effects on thermoelectronic transport properties", PHYSICAL REVIEW B89.
Example
julia> m = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0]; # effectuve mass
julia> band = ParabBand(m,0.0,1,1); # build the band
julia> μ = collect(-.1:.01:.1); # fermi level position
julia> model = BandStructure(1,band_1,μ) # build the band structure
julia> T = collect(50:10:650); # temperature range
julia> τ_form = Scattering.constant() # relaxation time
julia> S = seebeck_coefficient(model,T,τ_form); # compute the Seebeck coefficient
Mstar2t.Transport.carrier_concentration
— Functioncarrier_concentration(bandstructure::BandStructure, Ts::Union{Vector{Float64},Float64,Vector{Int64},Int64}, τ::Union{ScModel,Matthiessen}; exportasdf::Bool=false)
Function that computes the carrier concentration for a given bandstructure
and a choice of the relaxation time $\tau$, for given values of chemical potential $\mu$ and temperature Ts
. The transport coefficient is returned as a matrix
of dimensions (length(T),length(μ))
. The boolean variable exportasdf
allows to return the calculations as a DataFrame
with all the parameters included.
Example
julia> m = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0]; # effectuve mass
julia> band = ParabBand(m,0.0,1,1); # build the band
julia> μ = collect(-.1:.01:.1); # fermi level position
julia> model = BandStructure(1,band_1,μ) # build the band structure
julia> T = collect(50:10:650); # temperature range
julia> τ_form = Scattering.constant() # relaxation time
julia> n = carrier_concentration(model,T,τ_form); # compute the carrier concentration
Mstar2t.Transport.thermal_conductivity
— Functionthermal_conductivity(bandstructure::BandStructure, Ts::Union{Vector{Float64},Float64,Vector{Int64},Int64},τ::Union{ScModel,Matthiessen}; exportasdf::Bool=false, fulltensor::Bool=false)
Function that computes the thermal conductivity in units of $W/K$ for a given bandstructure
and a choice of the relaxation time $\tau$, for given values of chemical potential $\mu$ and temperature Ts
. The transport coefficient is returned as a matrix
of dimensions (length(T),length(μ))
. The boolean variable exportasdf
allows to return the calculations as a DataFrame
with all the parameters included. If fulltensor
is set to true
, the full tensor is returned in place of the trace.. Reference: "Theory of band warping and its effects on thermoelectronic transport properties", PHYSICAL REVIEW B89.
Example
julia> m = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0]; # effectuve mass
julia> band = ParabBand(m,0.0,1,1); # build the band
julia> μ = collect(-.1:.01:.1); # fermi level position
julia> model = BandStructure(1,band_1,μ) # build the band structure
julia> T = collect(50:10:650); # temperature range
julia> τ_form = Scattering.constant() # relaxation time
julia> Κ = thermal_conductivity(model,T,τ_form); # compute the thermal conductivity
Mstar2t.Transport.lorenz_tensor
— Functionlorenz_tensor(bandstructure::BandStructure, Ts::Union{Vector{Float64},Float64,Vector{Int64},Int64},τ::Union{ScModel,Matthiessen}; exportasdf::Bool=false, fulltensor::Bool=false)
Function that computes the Lorentz tensor for a given bandstructure
and a choice of the relaxation time $\tau$, for given values of chemical potential $\mu$ and temperature Ts
. The transport coefficient is returned as a matrix
of dimensions (length(T),length(μ))
. The boolean variable exportasdf
allows to return the calculations as a DataFrame
with all the parameters included. If fulltensor
is set to true
, the full tensor is returned in place of the trace. The Lorentz tensor is defined as the ratio between thermal and electrical conductivity multiplied by temperature. Ref.: https://en.wikipedia.org/wiki/Wiedemann-Franz_law
Example
julia> m = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0]; # effectuve mass
julia> band = ParabBand(m,0.0,1,1); # build the band
julia> μ = collect(-.1:.01:.1); # fermi level position
julia> model = BandStructure(1,band_1,μ) # build the band structure
julia> T = collect(50:10:650); # temperature range
julia> τ_form = Scattering.constant() # relaxation time
julia> L = lorenz_tensor(model,T,τ_form); # compute the lorenz tensor
Mstar2t.Transport.Ln
— FunctionLn(ϵ₀::Float64, μ::Float64, bandtype::Int64, T::Float64, β::Float64, n::Int64, idx::Int64, scm::ScModel)
This function is called when an integration is required. $n=0,1,2$ is the index of the kinetic coefficient. Reference: "Theory of band warping and its effects on thermoelectronic transport properties", PHYSICAL REVIEW B89.
Mstar2t.Scattering.constant
— Functionconstant(A::T=1.0) where {T<:Real}
Model for the relaxation time is set to the constant relaxation time approximation (CRTA), $\tau = A\cdot 10^{-14}$. Tje default value for A is 1.0.
Parameter:
A
: value for the A constant above.
Example
julia> τ_form = Scattering.constant() # relaxation time
Mstar2t.Scattering.T_fun
— FunctionT_fun(τ::Function)
Set the temperature dependence of the relaxation time from a function given by the user.
Example
julia> f(T) = sqrt(T) # τ ∝ √T
julia> τ_form = T_fun(f) # relaxation time
Mstar2t.Scattering.ET_fun
— FunctionET_fun(τ::Function)
Set the energy and temperature dependence of the relaxation time from a function given by the user.
Example
julia> f(ϵ,T) = sqrt(ϵ)/T # τ ∝ √ϵ/T
julia> τ_form = ET_fun(f) # relaxation time
Mstar2t.Scattering.impurity
— Functionimpurity(ϵ_im::Real, A_im::Real=1;γ::Real=1.)
Include impurity scattering in the simulations.
Parameters:
ϵ_im
: energy of the impurity in eVA
: multiplicative constant in front of the functional expressionγ
: γ-parameter three-parameter Lorentzian function (Ref: https://en.wikipedia.org/wiki/Cauchy_distribution)
Example
julia> ϵ_im = 0.2 # eV
julia> τ_form = impurity(ϵ_im, A=1e-1) # relaxation time
Mstar2t.Scattering.acoustic
— Functionacoustic(bands_min::Union{BandStructure,Real}, A_sm::Real=1., τm_max::Real=1.; T₀::Real=50., μ_min::Real=2, μ_max::Real=2, s=2)
Include acoustic scattering in the simulations. Implementation of the functional expression in Wilson, Alan Herries. 1953, The theory of metals by A. H. Wilson Cambridge Uni. Press Cambridge, England.
Parameters:
bands_min
`: energy of the lowest band. If a BandStructure type is passed, the value is automatically derived.A_sm
: multiplicative constant in front of the functional expression for the acoustic τ in semiconductors in units of $5\cdot 10^{-20}$ s.τm_max
: free parameter to constraint the functional form for the acoustic τ in metals. $τ(ϵ₀+μ\_max,T) = τm\_max$ in units of $2\cdot 10^{-12}$ s, where $ϵ₀$ is the minimum/maximum of the band.T₀
: minimum of temperature range in which the τ is defined (note: τ ∝ 1/(T-T₀), where T is the temperature input defined by the user.μ_min
: left shift in energy from the energy of the lowest band (ϵ_min
). It defines the first point at which τ is computed (i.e., $τ \propto 1/\sqrt{μ-(ϵ\_min-μ\_min)}$).μ_max
: right shift in energy from the energy of the band $\ϵ₀$. It defines the last point at which τ is computed (i.e., $τ(ϵ₀+μ\_max,T) = τm\_max$).
Example
julia> band_1 = ParabBand([5.0, 5.0, 5.0, 0.0, 0.0, 0.0],1.0, 1,1); # conduction band
julia> band_2 = ParabBand([0.1, 0.5, 3.0, 0.0, 0.0, 0.0],0.5,-1,1); # valence band
julia> model = BandStructure(2,[band_1,band_2],0.8); # build the two-band structure
julia> τ_form = Scattering.acoustic(model,T₀=180,μ_min=5,μ_max=5); # relaxation time
Mstar2t.Scattering.matthiessen
— Functionmatthiessen(τ_models::Array{ScModel};γ::Float64=-1.)
This function applies Matthiessen's rule to sum up different scattering mechanisms.
Parameters:
τ_models
: vector of relaxation time modelsγ
: exponent in the generalized Matthiessen's rule:
\[\tau = (\sum_i τ_i^γ)^{1/γ}, \]
where each $\tau_i$ can be a function of $T$ and/or $\mu$.
Example
julia> ϵ_im = 0.2 # eV
julia> τ1 = constant() # CRTA
julia> τ2 = impurity(ϵ_im,A=1e-1) # impurity scattering
julia> τ_model = matthiessen([τ1,τ2]); # matthiessen's rule
Mstar2t.Plot.plot
— Functionplot(num_plots, x_axis, data, z; vlines=Float64[], annotations=[], color=ColorTypes.RGBA{Float64}(Colors.JULIA_LOGO_COLORS.blue), colorscheme=:viridis, titles=L"", xlabels=L"", ylabels=L"", zlabel=L"")
Plot the transport coefficients. y
must be a vector if plotting a single line, or a matrix of shape (length(x_axis),length(z_axis))
if plotting more than one line.
Parameters:
- num_plots::Int64 : number of transport coefficients to plot in the same figure
- x::Union{Vector{Float64},Float64,Vector{Int64},Int64} # x axis vector
- y::Union{Vector{Float64},Matrix{Float64},Vector{Vector{Float64}},Vector{Matrix{Float64}}} # y axis vector/matrix
- z::Union{Vector{Float64},Vector{Int64}} # z axis vector
- titles::Union{LaTeXString,Vector{LaTeXString}} # titles for each plot
- xlabels::Union{LaTeXString,Vector{LaTeXString}} # x labels for each plot
- ylabels::Union{LaTeXString,Vector{LaTeXString}} # x labels for each plot
- zlabel::LaTeXString # z labels for the shared colorbar
- color::ColorTypes.RGBA{Float64} # line color
- colorscheme # colormap
- vlines::Union{Float64,Array{Float64}} # add vertical lines to the plot
- annotations::Vector{Any} # add text annotation to the plot
Example
Single line
(z is empty). Plot of electrical conductivity, Seebeck and carrier concentration as functions of the bandgap
.
julia> titles = [L"$σ$ vs gap, $μ$ = %$μ, $τ = const$", L"$S$ vs gap, $μ$ = %$μ, $τ = const$", L"$n$ vs gap, $μ$ = %$μ, $τ = const$"]
julia> xlabels = [L"$Gap$ $[eV]$" for i in 1:3];
julia> ylabels = [L"$\sigma$ $[(\Omega m)^{-1}]$", L"$S$ $[\mu VK^{-1}]$", L"n"];
julia> annotations = Dict(L"v1" => [(0.11,6e5),(0.11,0),(0.11,5e6)], L"v2" => [(0.01,6e5-5e3),(0.01,0),(0.01,5e6-5e4)]);
julia> fig = plot(3, gap, [σ,S,n], titles=titles, xlabels=xlabels, ylabels=ylabels, vlines=[0.0,0.1], annotations=annotations);
Multi-lines
. Plot of electrical conductivity, Seebeck and carrier concentration as functions of the band
gapand
temperature` (z axis).
julia> titles = [L"$σ$ vs gap, $μ$ = %$μ, $τ = const$", L"$S$ vs gap, $μ$ = %$μ, $τ = const$", L"$n$ vs gap, $μ$ = %$μ, $τ = const$"]
julia> xlabels = [L"$Gap$ $[eV]$" for i in 1:3];
julia> ylabels = [L"$\sigma$ $[(\Omega m)^{-1}]$", L"$S$ $[\mu VK^{-1}]$", L"n"];
julia> zlabel = L"$T$ $[K]$";
julia> annotations = Dict(L"v1" => [(0.11,6e5),(0.11,0),(0.11,5e6)], L"v2" => [(0.01,6e5-5e3),(0.01,0),(0.01,5e6-5e4)]);
julia> fig = plot(3, gap, [σ,S,n], T, titles=titles, xlabels=xlabels, ylabels=ylabels, zlabel=zlabel, vlines=[0.0,0.1], annotations=annotations);
`
Mstar2t.Plot.plot_bandstructure
— Functionplot_bandstructure(bs::BandStructure,xaxis::AbstractArray=range(-1, 1, length=100); colors=nothing, label="")
Plot a given parabolic band structure.
Parameters:
bs
: band structurexaxis
: x axis vectorcolors
: tuple of three color fromColors
. First: conduction bands. Second: valence band. Third: Fermi level.
Examples:
julia> band_1 = ParabBand([5.0, 5.0, 5.0, 0.0, 0.0, 0.0],1.0, 1,1); # conduction band
julia> band_2 = ParabBand([0.1, 0.5, 3.0, 0.0, 0.0, 0.0],0.5,-1,1); # valence band
julia> μ = 0.8; # Fermi level
julia> model = BandStructure(2,[band_1,band_2],μ); # build the two-band structure
julia> fig,ax = plot_bandstructure(model,colors=(:green,:red,:blue))
Mstar2t.Plot.plot_τ
— Functionplot_τ(scm::Union{ScModel,Matthiessen}, t::Union{Vector{Float64},Vector{Int64}}, type="e-T",E_argmax::Int64=50; μ=nothing, ϵ₀::Float64=-42., bandtype::Int64=0)
3D plot of a given relaxation time as a function of temperature, energy and Fermi level.
Parameters:
scm
: relaxation time modelt
: temperaturetype
: type of plot specificed in "xaxis-zaxis" format. Available types are "e-T","T-e","μ-T","T".μ
: Fermi level vector. It's mandatory if type="μ-T".ϵ₀
: band energy. It's mandatory for acoustic relaxation time.bandtype
: band type (conduction, valence). It's mandatory for acoustic relaxation time.
Examples:
- Constant relaxation time
julia> T = collect(300:10:650); # temperature
julia> τ_form = Scattering.constant()
julia> fig = plot_τ(τ_form, T, "T");
- Acoustic
julia> ϵ₀ = .1; # band energy
julia> type = 1; # band type (conduction)
julia> T = collect(300:10:650); # temperature
julia> μ = collect(0.0:0.1:1.0); # Fermi level
julia> τ_form = Scattering.acoustic(model,T₀=180,μ_min=5,μ_max=5);
julia> fig = plot_τ(τ_form, T, "μ-T", μ=μ, ϵ₀=ϵ₀, bandtype=type);
Mstar2t.Plot.savefig
— Functionsavefig(fullpath::String, fig::Figure)
Export to disk a given figure to fullpath
.
Example
- Save the plot for the relaxation time as a function of temperature and Fermi level.
julia> τ_form = Scattering.acoustic(model,T₀=180,μ_min=5,μ_max=5);
julia> fig = plot_τ(τ_form, T, "μ-T", μ=μ, ϵ₀=ϵ₀, bandtype=type);
julia> savefig("acoustic_tau.png", fig);
- Save the plot for the electrical conductivity, Seebeck coefficient and carrier density as a function of temperature and Fermi level.
julia> fig = plot(3, T, [σ,S*1e6,n], μ, titles=titles, xlabels=xlabels, ylabels=ylabels, zlabel=zlabel);
julia> savefig("transport_coefficients.png", fig);
Mstar2t.Utils.savedata
— Functionsavedata(path::String, data::DataFrame)
Export to disk a given tranport coefficients simulation to a csv file.