QCLayers module

This module contains QCLayers and SchrodingerLayer classes. The physics model can be seen in Quantum Tab

SchrodingerLayer is a Schrodinger solver together with transition properties evaluation, QCLayers is a wrapper of SchrodingerLayer with material information coded.

QCLayers class

class ErwinJr2.QCLayers.QCLayers(substrate='InP', materials=['InGaAs', 'AlInAs'], moleFracs=[0.53, 0.52], xres=0.5, Eres=0.5, statePerRepeat=20, layerWidths=[10.0], layerMtrls=None, layerDopings=None, customIFR=False, mtrlIFRLambda=None, mtrlIFRDelta=None, ifrDelta=None, ifrLambda=None, layerARs=None, EField=0, repeats=3, T=300.0, solver='ODE', description='', wl=3.0)

Class for Quantum Cascade Layers

Parameters
  • substrate (str) –

    The substrate material for the device, which determines the well and barrier material

    substrate

    well

    barrier

    InP

    InxGa1-xAs

    Al1-xInxAs

    GaAs

    AlxGa1-xAs

    AlxGa1-xAs

    GaSb

    InAsySb1-y

    AlxGa1-xSb

  • materials (List[str]) – Name of alloys for the heterostructure materials, len >= 2

  • moleFracs (List[float]) – mole fraction for each possible layer material, len = Mp. of materials

  • xres (float) – Position resolution, in Armstrong

  • Eres (float) – Energy resolution, in meV. This number being too large may results in miss of some states while this being too small will make a long computation time. The parameter does not mean the accuracy of the eigen energy. It’s required for algorithm reasons because of lack of a universal global root finding.

  • statePerRepeat (int) – Number of states per repeat, used for calculating matrixEigenCount

  • wl (float) – The wavelength for the design, in unit um, gain spectrum and optimization, but doesn’t go into quantum solver

  • layerWidths (List[float]) – Width of each layer, in angstrom. len = No. of layers

  • layerMtrls (List[int]) – Label of materials, depending on substrate. len = No. of layers

  • layerDopings (List[float]) – Doping per volume in unit 1e17 cm-3. len = No. of layers

  • layerARs (List[bool]) – Binaries indicating if the layer is active(True) or not(False), only affects basis solver. len = No. of layers

  • customIFR (bool) – Wether to use a customized IFR parameter rather than a material determined parameter.

  • mtrlIFRLambda (List[float]) – The interface roughness lambda after materials[n], len = No. of materials

  • mtrlIFRDelta (List[float]) – The interface roughness delta after materials[n], len = No. of materials

  • EField (float) – External (static) electrical field, in kV/cm = 1e5 V/m

  • repeats (int) – Number of repeat times for the given structure

  • T (float) – Temperature of the device, affecting material property

  • basisAROnly (bool) – For basis solver if only the Active Region (AR) should be solved.

  • basisARInjector (bool) – For basis solver if there should be separator between AR->Injector

  • basisInjectorAR (bool) – For basis solver if there should be separator between Injector->AR

  • solver (str) – The solver used for the eigen problem: ‘ODE’ or ‘matrix’. By default ‘ODE’ if C library exists, ‘matrix’ is a full back.

  • includeIFR (bool) – Weather to include IFR scattering for performance estimation.

  • matrixEigenCount (int) – The number of eigen pairs to calculate in the ‘matrix’ solver. It would be very expensive to calculate all of them.

  • status (str) –

    • ‘unsolved’ meaning the structure is not solved yet.

    • ’basis’ meaning the eigen problem is solved for basis

    • ’solved’ meaning the eigen problem is solved.

    • ’solved-full’ meaning the population distribution is known.

  • description (str) – Description of the data. For book-keeping purposes.

add_mtrl(mtrl=None, moleFrac=None, IFRLambda=None, IFRDelta=None)

Add a new material possibility

del_mtrl(n)

Delete materials labeled n. All layers of this material will become previous materials[n-1 if n >0 else 1]. There should be at least two materials otherwise there will be error.

effective_ridx(wl)

Return the effective refractive index for TM mode

Return type

Union[float, ndarray]

figure_of_merit(upper, lower)

Calculate Figure of Merit. This function must be called after solving for wave functions

Parameters
  • upper (int) – define the transition from upper to lower

  • lower (int) – define the transition from upper to lower

Return type

float

Returns

float (Figure of Merit)

Yields
  • tauLO_l (float) – the lower state lifetime from LO scattering

  • tauLO_u (float) – the upper state lifetime from LO scattering

  • tauLO_ul (float) – the transition rate from upper to lower due to LO scattering

  • tauIFR_l (float) – the lower state lifetime from IFR scattering

  • tauIFR_u (float) – the upper state lifetime from IFR scattering

  • tauIFR_ul (float) – the transition rate from upper to lower due to IFR scattering

  • tau_u (float) – 1/(1/tauLO_u + 1/tauIFR_u)

  • tau_l (float) – 1/(1/tauLO_l + 1/tauIFR_l)

  • FoM (float) – the Figure of Merit in unit angstrom^2 ps

full_gain_spectrum(wl=None)

Perform fully automatic calculation for the gain on wavelength(s).

Return type

Union[float, ndarray]

full_population()

Apart from SchrodingerLayer.full_population, this also yield current density, with knowledge of doping

Return type

array

layerMc(n)

The conduction band effective mass at n-th layer, in m0

Return type

float

layerVc(n)

The conduction band offset at n-th layer in eV

Return type

float

property mtrlOffset: float

Return the conduction band offset (difference between highest conduction band and lowest conduction band energy) of materials, in unit eV

Return type

float

property netStrain: float

Return average strain perpendicular to the layer plane, in percentage.

Return type

float

populate_material()

Update following band structure parameters (with type np.array of float): xVc, xVX, xVL, xVLH, xVSO, bandParams = (xEg, xF, xEp, xESO)

set_mtrl(n, mtrl=None, moleFrac=None)

Set material[n] to new material (mtrl) and/or moleFrac

property sheet_density: float

Return the sheet density of doping per period, in unit cm^-2

Return type

float

update_mtrls()

Update properties for the materials. This should be called every time the material parameters are changed directly via materials, moleFracs or temperature.

Yields
  • a_parallel (float) – The parallel crystal constant of the structure, determined by the substrate material.

  • mtrlAlloys (List[Material.Alloy]) – list of the Alloy objects for processing material properties.

SchrodingerLayer class

class ErwinJr2.QCLayers.SchrodingerLayer(xres=0.5, Eres=0.5, statePerRepeat=20, layerWidths=[10.0], layerARs=None, layerVc=None, layerMc=None, ifrDelta=None, ifrLambda=None, avghwLO=0.035, epsrho=69.0, EField=0.0, repeats=1)

Class for layer structure for Schrodinger solver.

This is used as the base class of QCLayers for separation of material property and the solver.

Parameters
  • xres (float) – Position resolution, in Armstrong

  • Eres (float) – Energy resolution, in meV. This number being too large may results in miss of some states while this being too small will make a long computation time. The parameter does not mean the accuracy of the eigen energy. It’s required for algorithm reasons because of lack of a universal global root finding.

  • statePerRepeat (int) – Number of states per repeat, used for calculating matrixEigenCount

  • layerWidths (List[float]) – Width of each layer, in angstrom. len = No. of layers

  • layerVc (Optional[List[float]]) – The conduction band offset of each layer, in unit eV, len = No. of layers

  • layerMc (Optional[List[float]]) – The conduction band effective mass of each layer, in unit m0, the free space electron mass, len = No. of layers

  • layerARs (List[bool]) – Binaries indicating if the layer is active(True) or not(False), only affects basis solver, len = No. of layers

  • ifrDelta (List[float]) – The standard deviation of the interface roughness for the interface at layer n and layer n+1, in unit angstrom, len = No. of layers. Default zero

  • ifrLambda (List[float]) – The correlation length of the interface roughness for the interface at layer n and layer n+1, in unit angstrom, len = No. of layers. Default zero

  • avghwLO (float) – The average LO phonon energy in unit eV. This value is used for LO phonon scattering calculation.

  • epsrho (float) – The effective relative permittivity for LO phonon. 1/epsrho = 1/(epsilon for high frequency) - 1/(epsilon for static) The value is used for LO phonon scattering calculation.

  • EField (float) – External (static) electrical field, in kV/cm = 1e5 V/m

  • repeats (int) – Number of repeat times for the given structure

  • crystalType (str) – Default being “simple”, meaning a simple parabolic effective mass is used for the calculation. For setting other than “simple”, populate_material should be implemented.

  • basisAROnly (bool) – For basis solver if only the Active Region (AR) should be solved.

  • basisARInjector (bool) – For basis solver if there should be separator between AR->Injector

  • basisInjectorAR (bool) – For basis solver if there should be separator between Injector->AR

  • solver (str) – The solver used for the eigen problem: ‘ODE’ or ‘matrix’. By default ‘ODE’ if C library exists, ‘matrix’ is a full back.

  • includeIFR (bool) – Weather to include IFR scattering for performance estimation.

  • matrixEigenCount (int) – The number of eigen pairs to calculate in the ‘matrix’ solver. It would be very expensive to calculate all of them.

  • status (str) –

    • ‘unsolved’ meaning the structure is not solved yet.

    • ’basis’ meaning the eigen problem is solved for basis

    • ’solved’ meaning the eigen problem is solved.

    • ’solved-full’ meaning the population distribution is known.

dephasing(upper, lower)

Calculate the broadening gamma of transition between upper -> lower transition, return gamma in unit eV as in Lorentzian:

\[\mathcal L(\omega) = \frac{1}{\pi} \frac{\gamma}{\gamma^2 + (\omega - \omega_0)^2}\]

If IFR scattering is included the broadening is calculated dominantly from IFR broadening and finite lifetime of upper and lower states. Otherwise 0.1 is returned.

Return type

float

full_population()

Calculate the electron full population on states, assuming the result of solve_whole and periodRecognize is valid and no state has coupling with states two or more periods away.

Return type

ndarray

Returns

population (np.array of float, dim = len(periodIdx)) – The population of electrons in the first recognized period, state label self.PeriodIdx[n]

Yields
  • transitions (np.array of float, dim = len(periodIdx)*len(periodIdx)) – The transition rate from self.PeriodIdx[i] to self.PeriodIdx[j]

  • decayRates (np.array of float, dim = len(periodIdx)) – Inverse of lifetimes

  • flow (float) – The flow of carrier, in unit ps^-1 (carrier density normalize to 1)

ifr_broadening(upper, lower)

Interface roughness induced broadening

Return type

float

ifr_lifetime(state)

Return to total life time due to IFR scattering. This is using cached results. if status is ‘solved-full’ and the state is a recognized state of a period, it’s calculated via translation of the wavefunction, otherwise it’s calculated based on row self.psis

Return type

float

ifr_transition(upper, lower)

Calculate the interface roughness (IFR) transition rate from upper to lower state at zero temperature, in unit ps^-1.

\[\frac{1}{\tau_{ij}^\text{IFR}} = \frac{\pi m^*_j}{\hbar^3} \sum_n \Delta_n^2\Lambda_n^2\delta U_n^2 \left|\psi_i(z_n)\psi_j^*(z_n)\right|^2 \mathrm e^{- \Lambda^2 m_j^* (E_i - E_j))/2\hbar^2}\]

This is using cached results. if status is ‘solved-full’ and the state is a recognized state of a period, it’s calculated via translation of the wavefunction, otherwise it’s calculated based on row self.psis

Return type

float

layerMc(n)

The conduction band effective mass at n-th layer, in m0

Return type

float

layerVc(n)

The conduction band offset at n-th layer in eV

Return type

float

lifetime(state)

A convenience wrap of return total lifetime of LO and IFR scattering or only LO scattering depending on self.includeIFR.

Return type

float

lo_lifetime(state)

Return the life time due to LO phonon scattering of the given state(label) This is using cached results. if status is ‘solved-full’ and the state is a recognized state of a period, it’s calculated via translation of the wavefunction, otherwise it’s calculated based on row self.psis

Return type

float

lo_transition(upper, lower)

The LO phonon transition lifetime from upper to lower, at zero temperature. This is using cached results. if status is ‘solved-full’ and the state is a recognized state of a period, it’s calculated via translation of the wavefunction, otherwise it’s calculated based on row self.psis

Return type

float

period_map_build(tol=0.03, etol=0.001)

Map states to self.singlePeriodIdx, self.periodMap[n] is a tuple of (state index in self.singlePeriodIdx, shift of period(s) or None meaning it’s not mapped.

Return type

List[Tuple[int, int]]

period_recognize(tol=5e-05)

Pick a set of eigen states as states in a period.

Return type

ndarray

Returns

singlePeriodIdx (np.array of int) – These are indices for the recognized states of a single period.

Yields

unBound (set of int) – includes index of states that are not well bounded.

populate_Kane_matrix()

Populate the finite difference Hamiltonian operation for the Kane 3 band model.

Return type

spmatrix

Returns

Hsparse (scipy.sparse.spmatrix) – A sparse matrix for the Hamiltonian

Yields
  • Hsparse (scipy.sparse.spmatrix) – Same above

  • Hbanded (np.ndarray) – Upper banded form of the matrix for the banded solver, as used in scipy.linalg.eig_banded

populate_material()

This should be overridden to yield bandParams

populate_x()

Calculate the properties in terms of position

Yields
  • xPoints (np.ndarray of float) – position grid

  • xLayerNums (np.ndarray of int) – at xPoints[i] it’s xLayerNums[i]-th layer

  • xVc (np.ndarray of float) – The band offset energy at each point

  • xMc (np.ndarray of float) – The effective mass at different points

psi_overlap(upper, lower, shift=0)

Return psi[upper] * psi[lower] with psi[lower] shifted by shift number of periods.

Return type

ndarray

reset_IFR_cache()

Reset cache for IFR scattering. This is exposed s.t. eigen solver results can be maintained after changing IFR setting.

solve_basis()

solve basis for the QC device, with each basis being the eigen mode of a separate part of the layer structure

Yields
  • eigenEs (np.array of float) – the eigenenergy of the layer structure

  • psis (np.array of float) – the wave function

Return type

ndarray

solve_whole()

Solve the whole structure. Will choose from _solve_whole_ode or _solve_whole_matrix according to self.solver.

Yields
  • eigenEs (np.array of float) – the eigenenergy of the layer structure

  • psis (np.array of float) – the wave function

Return type

ndarray

state_population(state)

This method is only valid after fullPopulation has been called

Return type

float

xLayerMask(n)

Return the mask for the given layer number n. A left and right extra point is included for plotting purposes.

Return type

ndarray

Others

ErwinJr2.QCLayers.BASISPAD = 100

padding barrier for basis solver, unit Angstrom

ErwinJr2.QCLayers.INV_INF = 1e-20

for infinite small decay rate (ns-1)

ErwinJr2.QCLayers.QCMaterial = {'GaAs': ['AlGaAs'], 'GaSb': ['InAsSb', 'AlGaSb'], 'InP': ['InGaAs', 'AlInAs']}

supported substrate/material set

exception ErwinJr2.QCLayers.StateRecognizeError(expression='')

Bases: Exception

Raised when QClayers cannot recognize a period of states. This typically means the number of periods is too small, or for single period structure the barrier at the end of the structure is not long enough to bound a state.

expression -- input expression in which the error occurred.
ErwinJr2.QCLayers.auto_gain(qcl, wls=None)

Perform automatic gain calculation from newly loaded a qcl object.

This is equivalent to

qcl.populate_x()
qcl.solve_whole()
qcl.period_recognize()
qcl.full_population()
result = qcl.full_gain_spectrum(wls)
ErwinJr2.QCLayers.optimize_global(qcl, iter=50)

A global optimization using gradient descent.

ErwinJr2.QCLayers.optimize_layer(qcl, n, upper, lower, iter=50)

Optimize FoM*Lorentzian for n-th layer thickness, assuming the state index does not change. optimization is performed by searching on the position resolution steps.

Warning: This cannot specify a correct state if there are states index crossing.

TODO: Use period recognizer to improve the algorithm

Example

Here is an example of how to construct a QCLayers class and solve for eigen states and wave functions. A sample json file is can be found here.

1with open("path/to/file.json") as f:
2    qcl = SaveLoad.qclLoad(f)
3
4qcl.layerSelected = 3
5qcl.populate_x()
6qcl.solve_whole()
7qcl.dipole(19, 15)
8qcl.FoM(19, 15)