Overview – Constants, Types, Constructors and Methods
Constants, Types and Constructors
ChargeTransport.BarrierLoweringOff
ChargeTransport.BarrierLoweringOn
ChargeTransport.BarrierLoweringType
ChargeTransport.BoundaryModelType
ChargeTransport.BulkRecombination
ChargeTransport.CalculationType
ChargeTransport.Data
ChargeTransport.Data
ChargeTransport.DiffusionEnhanced
ChargeTransport.ExcessChemicalPotential
ChargeTransport.ExcessChemicalPotentialGraded
ChargeTransport.FluxApproximationType
ChargeTransport.GeneralizedSG
ChargeTransport.GenerationBeerLambert
ChargeTransport.GenerationModelType
ChargeTransport.GenerationNone
ChargeTransport.GenerationUniform
ChargeTransport.GenerationUserDefined
ChargeTransport.InEquilibrium
ChargeTransport.InterfaceModelType
ChargeTransport.InterfaceNone
ChargeTransport.InterfaceRecombination
ChargeTransport.MixedOhmicSchottkyContact
ChargeTransport.ModelType
ChargeTransport.OhmicContact
ChargeTransport.OhmicContactModelType
ChargeTransport.OutOfEquilibrium
ChargeTransport.OuterBoundaryModelType
ChargeTransport.Params
ChargeTransport.Params
ChargeTransport.ParamsNodal
ChargeTransport.ParamsNodal
ChargeTransport.QType
ChargeTransport.SRHWithTrapsType
ChargeTransport.SRHWithoutTrapsType
ChargeTransport.ScharfetterGummel
ChargeTransport.ScharfetterGummelGraded
ChargeTransport.SchottkyContact
ChargeTransport.StandardFuncSet
ChargeTransport.Stationary
ChargeTransport.System
ChargeTransport.System
ChargeTransport.Transient
ChargeTransport.TrapCarrier
Methods
ChargeTransport.Blakemore
ChargeTransport.Boltzmann
ChargeTransport.FermiDiracMinusOne
ChargeTransport.FermiDiracOneHalfBednarczyk
ChargeTransport.FermiDiracOneHalfTeSCA
ChargeTransport.RHSContinuityEquations!
ChargeTransport.RHSPoisson!
ChargeTransport.add_trap_density_Poisson!
ChargeTransport.bflux!
ChargeTransport.breaction!
ChargeTransport.breaction!
ChargeTransport.breaction!
ChargeTransport.breaction!
ChargeTransport.bstorage!
ChargeTransport.build_system
ChargeTransport.build_system
ChargeTransport.charge_density
ChargeTransport.charge_density
ChargeTransport.degenerateLimit
ChargeTransport.electroNeutralSolution
ChargeTransport.enable_ionic_carrier!
ChargeTransport.enable_trap_carrier!
ChargeTransport.etaFunction
ChargeTransport.etaFunction
ChargeTransport.etaFunction
ChargeTransport.etaFunction!
ChargeTransport.etaFunction!
ChargeTransport.flux!
ChargeTransport.generic_operator!
ChargeTransport.get_BEE!
ChargeTransport.get_DOS!
ChargeTransport.get_current_val
ChargeTransport.get_current_val
ChargeTransport.get_density
ChargeTransport.get_density
ChargeTransport.get_density!
ChargeTransport.get_density!
ChargeTransport.get_density!
ChargeTransport.plotDiffusionEnhancements
ChargeTransport.plotDistributions
ChargeTransport.plot_IV
ChargeTransport.plot_densities
ChargeTransport.plot_doping
ChargeTransport.plot_doping
ChargeTransport.plot_electroNeutralSolutionBoltzmann
ChargeTransport.plot_energies
ChargeTransport.plot_energies
ChargeTransport.plot_solution
ChargeTransport.printJacobi
ChargeTransport.reaction!
ChargeTransport.reaction!
ChargeTransport.reaction!
ChargeTransport.set_bulk_recombination
ChargeTransport.set_plotting_labels
ChargeTransport.storage!
ChargeTransport.storage!
ChargeTransport.trap_density!
ChargeTransport.zeroVoltage
Description of Constant, Types, Constructors and Methods
Constants, Types and Constructors
ChargeTransport.BarrierLoweringType
— TypePossible types for barrier lowering model.
ChargeTransport.BoundaryModelType
— TypePossible types of boundary models.
ChargeTransport.CalculationType
— TypePossible types for calculation type.
ChargeTransport.FluxApproximationType
— TypePossible types of flux discretization schemes.
ChargeTransport.GenerationModelType
— TypePossible types for generation model.
ChargeTransport.InterfaceModelType
— TypePossible Types of interface model (interior boundary conditions).
ChargeTransport.ModelType
— TypePossible types which indicate, if we consider stationary or transient problem.
ChargeTransport.OhmicContactModelType
— TypePossible mathematical types of ohmic contact boundary model.
ChargeTransport.OuterBoundaryModelType
— TypePossible types of outer boundary model.
ChargeTransport.QType
— TypeType of charge carriers and the electric potential (corresponding to VoronoiFVM.jl).
ChargeTransport.SRHWithTrapsType
— TypePossible types for SRH recombination without traps.
ChargeTransport.SRHWithoutTrapsType
— TypePossible type for SRH recombination without traps.
ChargeTransport.StandardFuncSet
— TypeType of statistics functions.
ChargeTransport.BarrierLoweringOff
— Typeabstract type BarrierLoweringOff
Abstract type for the neglection of Schottky barrier lowering as boundary conditions.
ChargeTransport.BarrierLoweringOn
— Typeabstract type BarrierLoweringOn
Abstract type for the choice of Schottky barrier lowering as boundary conditions.
ChargeTransport.BulkRecombination
— Typemutable struct BulkRecombination
A struct holding all necessary information for building bulk recombination. With help of this constructor we can read out the indices the user chooses for electron and hole quasi Fermi potentials.
iphin::Int64
: Index for FVM construction of electron quasi Fermi potential.
iphip::Int64
: Index for FVM construction of hole quasi Fermi potential.
bulk_recomb_Auger::Bool
: Boolean for present Auger recombination in bulk.
bulk_recomb_radiative::Bool
: Boolean for present radiative recombination in bulk.
bulk_recomb_SRH::Union{Type{SRHOff}, Type{ChargeTransport.SRHStationary}, Type{SRHTrapsStationary}, Type{SRHTrapsTransient}}
: DataType for present SRH recombination in bulk. This needs to be a Type due to cases with or without mobile traps.
SRH_2species_trap::Union{Type{SRH2SpeciesPresentTrapDens}, Type{SRHOff}, Type{ChargeTransport.SRHStationary}}
: Data type with which you can include a stationary trap density to the right-hand side of the Poisson equation. This stationary trap density corresponds to the number of unoccupied trap states.
ChargeTransport.Data
— Typemutable struct Data{TFuncs<:Function, TVoltageFunc<:Function, TGenerationData<:Union{Array{Float64, 3}, Function, VecOrMat{Float64}}}
A struct holding all data information including model and numerics information, but also all physical parameters for a drift-diffusion simulation of a semiconductor device.
F::Vector{TFuncs} where TFuncs<:Function
: An array with the corresponding distribution function $\mathcal{F}_\alpha$ for all carriers $\alpha$.
qFModel::Union{Type{ContQF}, Type{DiscontQF}}
: An datatype containing the information, whether at least on quasi Fermi potential is assumend to be continuous or discontinuous.
boundaryType::Vector{Union{Type{InterfaceNone}, Type{InterfaceRecombination}, Type{MixedOhmicSchottkyContact}, Type{OhmicContact}, Type{SchottkyBarrierLowering}, Type{SchottkyContact}}}
: An array of DataTypes with the type of boundary model for each boundary (interior and exterior).
contactVoltageFunction::Vector{TVoltageFunc} where TVoltageFunc<:Function
: An array containing predefined functions for the applied bias in dependance of time at each outer boundary.
bulkRecombination::BulkRecombination
: A struct containing information concerning the bulk recombination model.
generationData::Union{Array{Float64, 3}, Function, VecOrMat{Float64}}
: A function/Array containing the user-specific photogeneration rate. It can be a function which is specified in the user example or an array which is read in and calculatd with, e.g., an external software.
isContinuous::Vector{Bool}
: An array containing information on whether charge carriers are continuous or discontinuous. This is needed for building the AbstractQuantities which handle the indices of charge carriers on different regions.
chargeCarrierList::Vector{Union{Int64, VoronoiFVM.ContinuousQuantity{Int32}, VoronoiFVM.DiscontinuousQuantity{Int32}, VoronoiFVM.InterfaceQuantity{Int32}}}
: This list stores all charge carriers with the correct type needed for VoronoiFVM.
electricCarrierList::Vector{Int64}
: This list stores all electric carrier indices, i.e. the one of electrons and holes.
ionicCarrierList::Vector{ChargeTransport.IonicCarrier}
: This list contains all defined ionic carriers as a struct of Type IonicCarrier with all needed information on the ionic carriers (can be either ions or ion vacancies).
trapCarrierList::Vector{ChargeTransport.TrapCarrier}
: This list contains all defined trap carriers for the SRH recombination as a struct of Type TrapCarrier with all needed information on the trap carriers.
AuxTrapValues::ChargeTransport.AuxiliaryStationaryTrapValues
: A struct which contains auxiliary trap values for the stationary setting.
index_psi::Union{Int64, VoronoiFVM.ContinuousQuantity{Int32}, VoronoiFVM.DiscontinuousQuantity{Int32}, VoronoiFVM.InterfaceQuantity{Int32}}
: This variable stores the index of the electric potential. Based on the user choice we have with this new type the opportunity to simulate discontinuous unknowns.
barrierLoweringInfo::ChargeTransport.BarrierLoweringSpecies
: This is a struct containing all information necessary to simulate Schottky Barrier Lowering.
fluxApproximation::Vector{Union{Type{DiffusionEnhanced}, Type{DiffusionEnhancedModifiedDrift}, Type{ExcessChemicalPotential}, Type{ExcessChemicalPotentialGraded}, Type{GeneralizedSG}, Type{ScharfetterGummel}, Type{ScharfetterGummelGraded}}}
: A DataType for the flux discretization method.
calculationType::Union{Type{InEquilibrium}, Type{OutOfEquilibrium}}
: A DataType for equilibrium or out of equilibrium calculations.
modelType::Union{Type{Stationary}, Type{Transient}}
: A DataType for transient or stationary calculations.
generationModel::Union{Type{GenerationBeerLambert}, Type{GenerationNone}, Type{GenerationUniform}, Type{GenerationUserDefined}}
: A DataType for for generation model.
λ1::Float64
: An embedding parameter used to solve the nonlinear Poisson problem, where for λ1 = 0 the right hand-side is set to zero whereas for for λ1 = 1 we have a full space charge density.
λ2::Float64
: An embedding parameter for the generation rate.
λ3::Float64
: An embedding parameter for an electrochemical reaction.
ohmicContactModel::Union{Type{OhmicContactDirichlet}, Type{OhmicContactRobin}}
: Possibility to change the implementation of the ohmic contact boundary model for the electric potential (Dirichlet or Robin)
tempBEE1::Vector{Float64}
: Within this template informations concerning the band-edge energy of each carrier is stored locally which saves allocations. We have two of such templates due to the two point flux approximation schemes.
tempBEE2::Vector{Float64}
: See the description of tempBEE1.
tempDOS1::Vector{Float64}
: Within this template informations concerning the effective DOS of each carrier is stored locally which saves allocations. We have two of such templates due to the two point flux approximation schemes.
tempDOS2::Vector{Float64}
: See the desciption of tempDOS2.
params::Params
: A struct holding all region dependent parameter information. For more information see struct Params.
paramsnodal::ParamsNodal
: A struct holding all space dependent parameter information. For more information see struct ParamsNodal.
ChargeTransport.Data
— MethodData(
grid,
numberOfCarriers;
contactVoltageFunction,
generationData,
statfunctions
) -> Data{StandardFuncSet, T, Vector{Float64}} where T<:Function
Simplified constructor for Data which only takes the grid and the numberOfCarriers as argument. Here, all necessary information including the physical parameters, but also some numerical information are located.
ChargeTransport.DiffusionEnhanced
— Typeabstract type DiffusionEnhanced
Abstract type for diffusion enhanced flux discretization, check M. Bessemoulin-Chatard, “A finite volume scheme for convection–diffusion equations with nonlinear diffusion derived from the Scharfetter–Gummel scheme”, Numerische Mathematik, vol. 121, pp. 637–670, 2012.
ChargeTransport.ExcessChemicalPotential
— Typeabstract type ExcessChemicalPotential
Abstract type for excess chemical potential flux discretization, check Z. Yu, and R. Dutton, “SEDAN III – A one-dimensional device simulator”, http://www-tcad.stanford.edu/tcad/programs/sedan3.html, 1988.
ChargeTransport.ExcessChemicalPotentialGraded
— Typeabstract type ExcessChemicalPotentialGraded
Abstract type for excess chemical potential flux discretization for graded effective density of states and/or graded band-edge energies. This means, use this flux when at least one of these parameters is assumed to be space-dependent.
ChargeTransport.GeneralizedSG
— Typeabstract type GeneralizedSG
Abstract type for generalized Scharfetter-Gummel flux discretization. This flux approximation results in an implicit equation which needs to be solved and is exact for all Blakemore type statistics functions with abritary γ, check T. Koprucki and K. Gärtner. “Discretization scheme for drift-diffusion equations with strong diffusion enhancement”. In: 12th International Conference on Numerical Simulation of Optoelectronic Devices (NUSOD). 2012, pp. 103–104.
ChargeTransport.GenerationBeerLambert
— Typeabstract type GenerationBeerLambert
Abstract type for Beer-Lambert generation. Note that this type is implemented, but not well tested yet.
ChargeTransport.GenerationNone
— Typeabstract type GenerationNone
Abstract type for no generation model.
ChargeTransport.GenerationUniform
— Typeabstract type GenerationUniform
Abstract type for uniform generation.
ChargeTransport.GenerationUserDefined
— Typeabstract type GenerationUserDefined
Abstract type for user defined generation.
ChargeTransport.InEquilibrium
— Typeabstract type InEquilibrium
Abstract type for equilibrium calculations.
ChargeTransport.InterfaceNone
— Typeabstract type InterfaceNone
Abstract type for no interface model.
ChargeTransport.InterfaceRecombination
— Typeabstract type InterfaceRecombination
Abstract type for surface recombination mechanisms.
ChargeTransport.MixedOhmicSchottkyContact
— TypeAbstract type for a mixed Ohmic and Schottky boundary model, resulting in all Dirichlet type conditions for electrons, holes and the electric potential.
ChargeTransport.OhmicContact
— TypeAbstract type for ohmic contacts as outer boundary model.
ChargeTransport.OutOfEquilibrium
— Typeabstract type OutOfEquilibrium
Abstract type for out of equilibrium calculations.
ChargeTransport.Params
— Typemutable struct Params
A struct holding the physical region dependent parameters for a drift-diffusion simulation of a semiconductor device.
numberOfNodes::Int64
: Number of nodes used for the disretization of the domain $\mathbf{\Omega}$.
numberOfRegions::Int64
: Number of subregions $\mathbf{\Omega}_k$ within the domain $\mathbf{\Omega}$.
numberOfBoundaryRegions::Int64
: Number of boundary regions $(\partial \mathbf{\Omega})_k$ such that $\partial \mathbf{\Omega} = \cup_k (\partial \mathbf{\Omega})_k$. Note that here are inner and outer boundaries calculated.
numberOfCarriers::Int64
: Number of moving charge carriers.
invertedIllumination::Int64
: Parameter for the direction of illumination. If illumination is coming from the left, then set this value to 1. Otherwise, if the illumination comes from the right, set this value to -1.
temperature::Float64
: A given constant temperature.
UT::Float64
: The thermal voltage, which reads $U_T = k_B T / q$.
γ::Float64
: The parameter of the Blakemore statistics (needed for the generalizedSG flux).
r0::Float64
: Prefactor of electro-chemical reaction of internal boundary conditions.
prefactor_SRH::Float64
: Prefactor for stationary SRH recombination.
generationPeak::Float64
: Parameter for the shift of generation peak of the Beer-Lambert generation profile.
SchottkyBarrier::Vector{Float64}
: An array for the given Schottky barriers at present Schotkky contacts.
contactVoltage::Vector{Float64}
: An array containing a constant value for the applied voltage.
bψEQ::Vector{Float64}
: An array containing a constant value for the electric potential in case of Dirichlet boundary conditions.
chargeNumbers::Vector{Float64}
: An array with the corresponding charge numbers $z_\alpha$ for all carriers $\alpha$.
bBandEdgeEnergy::Matrix{Float64}
: An array with the corresponding boundary band-edge energy values $E_\alpha$ in each region for each carrier $\alpha$.
bDensityOfStates::Matrix{Float64}
: An array with the corresponding boundary effective density of states values $N_\alpha$ for each carrier $\alpha$.
bMobility::Matrix{Float64}
: A 2D array with the corresponding boundary mobility values $\mu_\alpha$ in each boundary region for each carrier $\alpha$.
bDoping::Matrix{Float64}
: A 2D array with the corresponding boundary doping values for each carrier $\alpha$.
bVelocity::Matrix{Float64}
: A 2D array with the corresponding boundary velocity values for each carrier $\alpha$, when assuming Schottky contacts.
bReactionCoefficient::Matrix{Float64}
: An array to define the reaction coefficient at internal boundaries.
recombinationSRHvelocity::Matrix{Float64}
: A 2D array with the corresponding recombination surface boundary velocity values for electrons and holes.
bRecombinationSRHTrapDensity::Matrix{Float64}
: A 2D array with the corresponding recombination surface boundary density values for electrons and holes.
bRecombinationSRHLifetime::Matrix{Float64}
: A 2D array with the corresponding recombination surface recombination velocities.
bDensityEQ::Matrix{Float64}
: A 2D array containing the equilibrium density of electric charge carriers at the boundary.
doping::Matrix{Float64}
: A 2D array with the corresponding doping values for each carrier $\alpha$ on each region.
densityOfStates::Matrix{Float64}
: A 2D array with the corresponding effective density of states values $N_\alpha$ for each carrier $\alpha$ on each region.
bandEdgeEnergy::Matrix{Float64}
: A 2D array with the corresponding band-edge energy values $E_\alpha$ for each carrier $\alpha$ on each region.
mobility::Matrix{Float64}
: A 2D array with the corresponding mobility values $\mu_\alpha$ for each carrier $\alpha$ on each region.
recombinationSRHLifetime::Matrix{Float64}
: A 2D array with the corresponding SRH lifetimes $\tau_n, \tau_p$ for electrons and holes.
recombinationSRHTrapDensity::Matrix{Float64}
: A 2D array with the corresponding time-independent SRH trap densities $n_{\tau}, p_{\tau}$ for electrons and holes.
recombinationAuger::Matrix{Float64}
: A 2D array with the corresponding Auger coefficients for electrons and holes.
dielectricConstant::Vector{Float64}
: A region dependent dielectric constant.
dielectricConstantImageForce::Vector{Float64}
: A region dependent image force dielectric constant.
generationIncidentPhotonFlux::Vector{Float64}
: A region dependent array for the prefactor in the generation process which is the incident photon flux.
generationUniform::Vector{Float64}
: A region dependent array for an uniform generation rate.
generationAbsorption::Vector{Float64}
: A region dependent array for the absorption coefficient in the generation process.
recombinationRadiative::Vector{Float64}
: A region dependent array for the radiative recombination rate.
ChargeTransport.Params
— MethodParams(grid, numberOfCarriers) -> Params
Simplified constructor for Params which only takes the grid and the numberOfCarriers as argument.
ChargeTransport.ParamsNodal
— Typemutable struct ParamsNodal
A struct holding the physical nodal, i.e. space-dependent parameters for a drift-diffusion simulation of a semiconductor device.
dielectricConstant::Vector{Float64}
: A node dependent dielectric constant.
doping::Vector{Float64}
: A 1D array with the corresponding doping values on each node.
mobility::Matrix{Float64}
: A 2D array with the corresponding mobility values $\mu_\alpha$ for each carrier $\alpha$ on each node.
densityOfStates::Matrix{Float64}
: A 2D array with the corresponding effective density of states values $N_\alpha$ for each carrier $\alpha$ on each node.
bandEdgeEnergy::Matrix{Float64}
: A 2D array with the corresponding band-edge energy values $E_\alpha$ for each carrier $\alpha$ on each node.
ChargeTransport.ParamsNodal
— MethodParamsNodal(grid, numberOfCarriers) -> ParamsNodal
Simplified constructor for ParamsNodal which only takes the grid and the numberOfCarriers as argument.
ChargeTransport.ScharfetterGummel
— Typeabstract type ScharfetterGummel
Abstract type for Scharfetter-Gummel flux discretization. Choose this one, when the Boltzmann statistics function is chosen as statistics, check D. Scharfetter and H. Gummel, “Large-signal analysis of a silicon Read diode oscillator”, IEEE Trans. Electr. Dev., vol. 16, pp. 64–77, 1969.
ChargeTransport.ScharfetterGummelGraded
— Typeabstract type ScharfetterGummelGraded
Abstract type for Scharfetter-Gummel flux discretization for graded effective density of states and/or graded band-edge energies. This means, use this flux when at least one of these parameters is assumed to be space-dependent.
ChargeTransport.SchottkyContact
— TypeAbstract type for schottky contacts as boundary model.
ChargeTransport.Stationary
— Typeabstract type Stationary
Abstract type for stationary simulations.
ChargeTransport.System
— Typemutable struct System
A struct holding all information necessary for a drift-diffusion type system.
data::Data
: A struct holding all data information, see Data
fvmsys::VoronoiFVM.AbstractSystem
: A struct holding system information for the finite volume system.
ChargeTransport.System
— MethodSystem(grid, data; kwargs...)
System constructor which builds all necessary information needed based on the input parameters with special regard to the quasi Fermi potential model. This is the main struct in which all information on the input data, but also on the solving system, are stored.
ChargeTransport.Transient
— Typeabstract type Transient
Abstract type for transient simulations.
ChargeTransport.TrapCarrier
— Typemutable struct TrapCarrier
A struct holding all information necessary for enabling traps in the SRH recombination. With help of this constructor we can read out the index the user chooses for trap quasi Fermi potentials and the respective regions in which they are defined.
trapCarrier::Int64
: Index of trap carrier user defines.
regions::Vector{Int64}
: Corresponding regions where trap carrier is assumed to be present.
Methods
ChargeTransport.Blakemore
— MethodBlakemore(x::Real, γ::Real) -> Any
The Blakemore approximation $1/(\exp(-x) + γ)$ with variable real scalar $γ$, see J. S. Blakemore. “The Parameters of Partially Degenerate Semiconductors”. In: Proceedings of the Physical Society. Section A 65 (1952), pp. 460–461.
ChargeTransport.Boltzmann
— MethodBoltzmann(x::Real) -> Any
The Boltzmann statistics function $\exp(x)$.
ChargeTransport.FermiDiracMinusOne
— MethodFermiDiracMinusOne(x::Real) -> Any
The Fermi-Dirac integral of order $-1$ which reads $1/(\exp(-x) + 1)$, see J.S. Blakemore, Approximations for Fermi-Dirac integrals, especially the function $F_{1/2} (\eta)$ used to describe electron density in a semiconductor, Solid-State Electronics 25 (11) (1982) 1067 – 1076.
ChargeTransport.FermiDiracOneHalfBednarczyk
— MethodFermiDiracOneHalfBednarczyk(x::Real) -> Any
The incomplete Fermi-Dirac integral of order 1/2, implemented according to [Bednarczyk1978, "The Approximation of the Fermi-Dirac integral $F_{1/2}(\eta)$"].
ChargeTransport.FermiDiracOneHalfTeSCA
— MethodFermiDiracOneHalfTeSCA(x::Real) -> Any
The incomplete Fermi-Dirac integral of order 1/2, implemented according to the software package TeSCA, see https://wias-berlin.de/software/index.jsp?lang=1&id=TeSCA.
Modified to use log1p(x)=log(1+x).
ChargeTransport.RHSContinuityEquations!
— MethodRHSContinuityEquations!(f, u, node, data)
Function which builds right-hand side of electric charge carriers.
ChargeTransport.RHSPoisson!
— MethodRHSPoisson!(f, u, node, data, ipsi) -> Any
Function which builds right-hand side of Poisson equation, i.e. which builds the space charge density.
ChargeTransport.add_trap_density_Poisson!
— Methodadd_trap_density_Poisson!(; data, zt, Nt)
This method includes traps for a simplified model, where the trap carriers are not considered as additional carrier with an own continuity equation. In this case the trap density is additionally added to the right-hand side of Poisson equation.
ChargeTransport.bflux!
— Methodbflux!(f, u, bedge, data)
Master bflux! function. This is the function which enters VoronoiFVM and hands over for each boundary the flux within the boundary.
ChargeTransport.breaction!
— Methodbreaction!(f, u, bnode, data) -> Any
Master breaction! function. This is the function which enters VoronoiFVM and hands over for each boundary the chosen boundary model.
ChargeTransport.breaction!
— Methodbreaction!(
f,
u,
bnode,
data,
_::Type{OhmicContactDirichlet}
)
Creates ohmic boundary conditions via Dirichlet BC for the electrostatic potential $\psi$
$\psi = \psi_0 + U$,
where $\psi_0$ contains some given value and $U$ is an applied voltage.
$f[\psi] = -q/\delta \sum_\alpha{ z_\alpha (n_\alpha - C_\alpha) },$
where $C_\alpha$ corresponds to some doping w.r.t. the species $\alpha$.
The boundary conditions for electrons and holes are dirichlet conditions, where
$\varphi_{\alpha} = U.$`
ChargeTransport.breaction!
— Methodbreaction!(f, u, bnode, data, _::Type{OhmicContactRobin})
Creates ohmic boundary conditions via a penalty approach with penalty parameter $\delta$. For example, the right-hand side for the electrostatic potential $\psi$ is implemented as
$f[\psi] = -q/\delta ( (p - N_a) - (n - N_d) )$,
assuming a bipolar semiconductor. In general, we have for some given charge number $z_\alpha$
$f[\psi] = -q/\delta \sum_\alpha{ z_\alpha (n_\alpha - C_\alpha) },$
where $C_\alpha$ corresponds to some doping w.r.t. the species $\alpha$.
The boundary conditions for electrons and holes are dirichlet conditions, where
$\varphi_{\alpha} = U$`
with $U$ as an applied voltage.
ChargeTransport.breaction!
— Methodbreaction!(
f,
u,
bnode,
data,
_::Type{SchottkyBarrierLowering}
) -> Any
Creates Schottky boundary conditions with additional lowering which are modelled as
$\psi = - \phi_S/q + \sqrt{ -\frac{ q \nabla_{\boldsymbol{\nu}} \psi_\mathrm{R}}{4\pi \varepsilon_\mathrm{i}}} + U$,
where $\psi_\mathrm{R}$ denotes the electric potential with standard Schottky contacts and the same space charge density as $\psi$ and where $\varepsilon_\mathrm{i}}}$ corresponds to the image force dielectric constant.
To solve for this additional boundary conditions the projected gradient $\nabla_{\boldsymbol{\nu}} \psi_\mathrm{R}$ is stored within a boundary species and calculated in the method generic_operator!().
ChargeTransport.bstorage!
— Methodbstorage!(f, u, bnode, data) -> Any
Master bstorage! function. This is the function which enters VoronoiFVM and hands over for each boundary the time-dependent part of the chosen boundary model.
ChargeTransport.build_system
— Methodbuild_system(
grid,
data,
::Type{ContQF};
kwargs...
) -> System
The core of the system constructor. Here, the system for continuous quasi Fermi potentials is build.
ChargeTransport.build_system
— Methodbuild_system(
grid,
data,
::Type{DiscontQF};
kwargs...
) -> System
The core of the system constructor. Here, the system for discontinuous quasi Fermi potentials is build.
ChargeTransport.charge_density
— Methodcharge_density(
psi0,
phi,
UT,
EVector,
chargeNumbers,
dopingVector,
dosVector,
FVector
) -> Any
Compute the charge density, i.e. the right-hand side of Poisson's equation.
ChargeTransport.charge_density
— Methodcharge_density(ctsys, sol) -> Any
Compute the charge density for each region separately.
ChargeTransport.degenerateLimit
— MethoddegenerateLimit(x) -> Any
Degenerate limit of incomplete Fermi-Dirac integral of order 1/2.
ChargeTransport.electroNeutralSolution
— MethodelectroNeutralSolution(ctsys) -> Any
Compute the electro-neutral solution for the Boltzmann approximation. It is obtained by setting the left-hand side in the Poisson equation equal to zero and solving for $\psi$. The charge carriers may obey different statitics functions. Currently, this one is not well tested for the case of charge carriers beyond electrons and holes.
ChargeTransport.enable_ionic_carrier!
— Methodenable_ionic_carrier!(data; ionicCarrier, regions)
This method takes the user information concerning present ionic charge carriers, builds a struct of Type IonicCarrier and add this struct to the ionicCarrierList.
ChargeTransport.enable_trap_carrier!
— Methodenable_trap_carrier!(; data, trapCarrier, regions)
This method takes the user information concerning present trap charge carriers, builds a struct of Type TrapCarrier and add this struct to the trapCarrierList.
ChargeTransport.etaFunction!
— MethodetaFunction!(u, bnode::VoronoiFVM.BNode, data, icc) -> Any
The argument of the statistics function for boundary nodes.
ChargeTransport.etaFunction!
— MethodetaFunction!(u, node::VoronoiFVM.Node, data, icc) -> Any
The argument of the statistics function for interior nodes.
ChargeTransport.etaFunction
— MethodetaFunction(psi, phi, UT, E, z) -> Any
The argument of the statistics function for given $\varphi_\alpha$ and $\psi$
$z_\alpha / U_T ( (\varphi_\alpha - \psi) + E_\alpha / q ).$
The parameters $E_\alpha$ and $z_\alpha$ are given as vectors. This function may be used to compute the charge density, i.e. the right-hand side of the Poisson equation.
ChargeTransport.etaFunction
— MethodetaFunction(
u,
data,
node,
region,
icc,
in_region::Bool
) -> Any
The argument of the distribution function for floats.
ChargeTransport.etaFunction
— MethodetaFunction(
sol,
ireg::Int64,
ctsys,
icc::Union{Int64, VoronoiFVM.ContinuousQuantity{Int32}, VoronoiFVM.DiscontinuousQuantity{Int32}, VoronoiFVM.InterfaceQuantity{Int32}}
) -> Any
The argument of the statistics function for a given solution on a given interior region.
ChargeTransport.flux!
— Methodflux!(f, u, edge, data) -> Any
Master flux functions which enters VoronoiFVM. Flux discretization scheme is chosen in two steps. First, we need to see, if we are in or out of equilibrium. If, InEquilibrium, then no flux is passed. If outOfEquilibrium, we choose the flux approximation which the user chose for each charge carrier. For the displacement flux we use a finite difference approach.
ChargeTransport.generic_operator!
— Methodgeneric_operator!(f, u, fvmsys) -> Any
Generic operator to save the projected gradient of electric potential (for system with standard Schotty contacts). Note that this currently only working in one dimension!
ChargeTransport.get_BEE!
— Methodget_BEE!(
icc::Union{Int64, VoronoiFVM.ContinuousQuantity{Int32}, VoronoiFVM.DiscontinuousQuantity{Int32}, VoronoiFVM.InterfaceQuantity{Int32}},
node::VoronoiFVM.Node,
data
) -> Any
Defining locally the band-edge energy for interior nodes (analougesly for boundary nodes and edges).
ChargeTransport.get_DOS!
— Methodget_DOS!(
icc::Union{Int64, VoronoiFVM.ContinuousQuantity{Int32}, VoronoiFVM.DiscontinuousQuantity{Int32}, VoronoiFVM.InterfaceQuantity{Int32}},
node::VoronoiFVM.Node,
data
) -> Any
Defining locally the effective DOS for interior nodes (analogously for boundary nodes and edges).
ChargeTransport.get_current_val
— MethodCalculates current for time dependent problem.
ChargeTransport.get_current_val
— MethodCalculates current for stationary problem.
ChargeTransport.get_density!
— Methodget_density!(u, bnode::VoronoiFVM.BNode, data, icc) -> Any
For given potentials, compute corresponding densities for interior nodes.
ChargeTransport.get_density!
— Methodget_density!(
u,
edge::VoronoiFVM.Edge,
data,
icc
) -> Tuple{Any, Any}
For given potentials, compute corresponding densities for edges.
ChargeTransport.get_density!
— Methodget_density!(u, node::VoronoiFVM.Node, data, icc) -> Any
For given potentials, compute corresponding densities for interior nodes.
ChargeTransport.get_density
— Methodget_density(sol, data, icc, ireg; inode) -> Any
The densities for given potentials $\varphi_\alpha$ and $\psi$
ChargeTransport.get_density
— Methodget_density(
sol,
ireg::Int64,
ctsys,
icc::Union{Int64, VoronoiFVM.ContinuousQuantity{Int32}, VoronoiFVM.DiscontinuousQuantity{Int32}, VoronoiFVM.InterfaceQuantity{Int32}}
) -> Any
For given potentials, compute corresponding densities for given interior region corresponding to a homogeneous set of parameters.
ChargeTransport.plotDiffusionEnhancements
— MethodplotDiffusionEnhancements(; Plotter)
Plot diffusion enhancements.
ChargeTransport.plotDistributions
— MethodplotDistributions(; Plotter)
Plot different distribution integrals.
ChargeTransport.plot_IV
— Methodplot_IV(
Plotter,
biasValues,
IV,
title;
plotGridpoints
) -> Any
Method for showing the total current. One input parameter is the boolean plotGridpoints which makes it possible to plot markers, which indicate where the nodes are located.
ChargeTransport.plot_densities
— Methodplot_densities(
Plotter,
ctsys,
solution,
title,
label_density;
plotGridpoints
) -> Any
Plotting routine, where the charge carrier densities are depicted in dependence of space. The case of heterojunctions is tested, but yet multidimensional plottings are not included. One input parameter is the boolean plotGridpoints which makes it possible to plot markers, which indicate where the nodes are located.
ChargeTransport.plot_doping
— Methodplot_doping(Plotter, ctsys, label_density) -> Any
Possibility to plot the considered doping. This is especially useful for making sure that the interior and the boundary doping agree.
ChargeTransport.plot_doping
— MethodPlot doping for nodal dependent doping.
ChargeTransport.plot_electroNeutralSolutionBoltzmann
— Methodplot_electroNeutralSolutionBoltzmann(
Plotter,
grid,
psi0;
plotGridpoints
) -> Any
Plotting routine for depicting the electroneutral potential. One input parameter is the boolean plotGridpoints which makes it possible to plot markers, which indicate where the nodes are located.
ChargeTransport.plot_energies
— Methodplot_energies(
Plotter,
ctsys,
solution,
title,
label_energy;
plotGridpoints
) -> Any
With this method it is possible to plot the energies
$E_\alpha - q \psi \quad \text{w.r.t. space.}$
The case of heterojunctions is tested, but yet multidimensional plottings are not included.
One input parameter is the boolean plotGridpoints which makes it possible to plot markers, which indicate where the nodes are located.
ChargeTransport.plot_energies
— Methodplot_energies(Plotter, ctsys, label_BEE)
With this method it is possible to depict the band-edge energies $E_\alpha$. This can be useful for debugging when dealing with heterojunctions.
ChargeTransport.plot_solution
— Methodplot_solution(
Plotter,
ctsys,
solution,
title,
label_solution;
plotGridpoints
) -> Any
Method for plotting the solution vectors: the electrostatic potential $\psi$ as well as the charge carriers. The case of heterojunctions is tested, but yet multidimensional plottings are not included. One input parameter is the boolean plotGridpoints which makes it possible to plot markers, which indicate where the nodes are located.
ChargeTransport.printJacobi
— MethodprintJacobi(node, sys)
First try of debugger. Print the Jacobi matrix for a given node, i.e. the number of node in the grid and not the excact coordinate. This is only done for the one dimensional case so far.
ChargeTransport.reaction!
— Methodreaction!(f, u, node, data)
Master reaction! function. This is the function which enters VoronoiFVM and hands over reaction terms for concrete calculation type and bulk recombination model.
ChargeTransport.reaction!
— Methodreaction!(f, u, node, data, _::Type{InEquilibrium})
Reaction in case of equilibrium, i.e. no generation and recombination is considered.
ChargeTransport.reaction!
— Methodreaction!(f, u, node, data, _::Type{OutOfEquilibrium})
Sets up the right-hand sides. Assuming a bipolar semiconductor the right-hand side for the electrostatic potential becomes $f[ψ] = - q ((p - N_a) - (n - N_d) ) = - q \sum n_\alpha (n_\alpha - C_\alpha)$ for some doping $C_\alpha$ w.r.t. to the species $\alpha$. The right-hand sides for the charge carriers read as $f[n_\alpha] = - z_\alpha q (G - R)$ for all charge carriers $n_\alpha$. The recombination includes radiative, Auger and Shockley-Read-Hall recombination. For latter recombination process the stationary simplification is implemented. The recombination is only implemented for electron and holes and assumes that the electron index is 1 and the hole index is 2.
ChargeTransport.set_bulk_recombination
— Methodset_bulk_recombination(
;
iphin,
iphip,
bulk_recomb_Auger,
bulk_recomb_radiative,
bulk_recomb_SRH
)
Corresponding constructor for the bulk recombination model.
ChargeTransport.set_plotting_labels
— Methodset_plotting_labels(
data
) -> Tuple{Any, Any, Matrix{String}, Any}
Method which can be used to construct the arrays parsed to the plotting routines for labeling. The description for electrons and holes are predefined. If one wishes to extend by labels for, e.g. mobile ionic carriers or traps, this can be done within the main file.
ChargeTransport.storage!
— Methodstorage!(f, u, node, data) -> Any
Master storage! function. This is the function which enters VoronoiFVM and hands over a storage term, if we consider transient problem.
ChargeTransport.storage!
— Methodstorage!(
f,
u,
node,
data,
_::Type{OutOfEquilibrium}
) -> Float64
The storage term for time-dependent problems. Currently, for the time-dependent current densities the implicit Euler scheme is used. Hence, we have $f[n_\alpha] = z_\alpha q ∂_t n_\alpha$ and for the electrostatic potential $f[ψ] = 0$.
ChargeTransport.trap_density!
— Methodtrap_density!(icc, ireg, params, Et)
Compute trap densities for a given trap energy. [Currently, only done for the Boltzmann statistics and for region dependent parameters.]
ChargeTransport.zeroVoltage
— MethodFunction in case of an applied voltage equal to zero at one boundary.