qtcad.device.device module

Device objects and methods.

class qtcad.device.device.Device(mesh=None, inp_dict=None, conf_carriers='unspecified')

Bases: device.core.device.Device

Device properties. All units are SI. Inherits the methods of the device class of the poissonfem package (2D).

Attributes
  • mesh (mesh object) – The mesh on which device properties are defined

  • boundaries (list) – List of boundary objects for the physical boundary conditions (Poisson’s equation).

  • bnd_schrod (list) – List of boundary objects for the physical boundary conditions (Schrodinger’s equation).

  • mc (1D field) – Density of states effective mass of electrons at each node. Rows: Elements Cols: Nodes

  • mv (2D field) – Density of states effective mass of holes at each node (1st index) and for each band (2nd index).

  • g_star (float) – Effective electron g-factor (for device). One g-tensor for entire device - not position dependent. To get the associated g-tensor use get_electron_g_factor

  • hole_Zeeman_params (1D array) – 2 entries kappa and q. One set of ‘g-factors’ used for the entire device - not position dependent. The first entry is the kappa g-factor and the second is the q g-factor. To get the associated g-tensors use get_hole_Zeeman_params

  • B (2D field) – Magnetic field at each node of device.

  • hole_kp_model (string) – Function of k.p parameters to use to calculate the hole D matrix. Available values are: ‘luttinger_kohn_foreman’ and None.

  • D (5D field or None) – Hole D matrix if holes are considered. Default is None.

  • eps (1D field) – Permittivity at each node.

  • T (float) – Temperature

  • EF (float) – Fermi energy. Absolute value.

  • ND (1D field) – Density of donor impurities at each node.

  • NA (1D field) – Density of acceptor impurities at each node.

  • gD (int) – Donor level degeneracy.

  • gA (int) – Acceptor level degeneracy.

  • gc (1D field) – Total (spin and valley) conduction band degeneracy

  • gv (2D field) – Total valence band degeneracy for each node (1st index) and band (2nd index).

  • gqc (1D field) – Total degeneracy (spin and valley) of electrons to use when populating energy eigenstates. Default: 4, for spin degeneracy and two-fold valley degeneracy of the Delta_2 valley in silicon.

  • chi (1D field) – Semiconductor electron affinity at each node.

  • Eg (1D field) – Semiconductor bandgap at each node.

  • Ed (1D field) – Donor binding energy at each node.

  • Ea (1D field) – Acceptor binding energy at each node.

  • phi_F (float or 1D field) – Reference potential.

  • V (1D field) – Confinement potential at each node.

  • Vext (1D field) – External (offset) potential at each node. Used for heterostructures or for testing. An external potential is an additional potential accounted for only in the Schrodinger solver.

  • insulator (1D field) – Boolean (True of False) for each node. If True, then no carrier can be present at the corresponding node.

  • dot_region (list) – Submeshes or string labels defining the quantum dot region, in which confined carrier densities are set to zero in non-linear Poisson

  • dot_nodes (1D array or None) – Index of each global node belonging to the quantum dot region, where confined carrier densities are set to zero in non-linear Poisson.

  • dot_elems (1D array or None) – Index of each element belonging to the quantum dot region, where confined carrier densities are set to zero in non-linear Poisson.

  • n (1D field) – Total electron density at each node, including contributions from both classical semiconductor physics and quantum effects.

  • p (1D field) – Total hole density at each node, including contributions from both classical semiconductor physics and quantum effects.

  • nq (1D field) – Electron density at each node due to quantum effects, i.e., the solution of Schrodinger’s equation.

  • pq (1D field) – Hole density at each node due to quantum effects, i.e., the solution of Schrodinger’s equation.

  • Nplus (1D field) – Density of ionized donors.

  • Nminus (1D field) – Density of ionized acceptors.

  • rho (1D field) – Total charge density at each node.

  • rho_prime (1D field) – Derivative of the charge density with respect to electrostatic potential phi at each node.

  • rho_0 (field) – Background charge density.

  • energies (1D array or None) – Eigen-energies of Schrodinger’s equation. For holes, the convention used is the one presented in QTCAD convention for holes.

  • eigenfunctions (2D field, 3D field or None) – Eigenfunctions of Schrodinger’s equation. For electrons, the first index is global nodes, the second is energy level. For holes, the first index is global nodes, the second is energy level, and the third is band. Moreover, for holes, the convention used is the one presented in QTCAD convention for holes.

  • macro_diff (1D field) – Difference between the macroscopic potential of the material and that of the reference material.

  • vlnce_band_macro (1D field) – Position of the valence band edge with respect to the macroscopic potential in the bulk material.

  • conf_carriers (str or None) – Which carrier densities from quantum-confined regions to add to the total charge density. Electrons: ‘electrons’ or ‘e’ Holes: ‘holes’ or ‘h’ Neither: None.

  • statistics (str) – Model to use for the carrier statistics. Maxwell-Boltzmann (‘MB’), Fermi-Dirac (‘FD’), or approximate Fermi-Dirac (‘FD_approx’). Default: ‘FD_approx’

  • ionization (str) – Dopant ionization model: ‘complete’ or ‘incomplete’.

  • coulomb_mat (4D array) – Coefficients \(V_{ijkl}\) of the Coulomb interaction in second quantization. Typically a 4D array, but can also be a 2D array when overlaps between basis states are neglected, in which case the array elements are related to the general 4D array by \(V_{ij} = V_{ijij}\).

  • many_body_subspaces (list of Subspace objects) – Subspace objects resulting from the many-body solver. The n-th element of the list contains information about the n-particle subspace.

  • chem_potentials (1D array) – Chemical potentials of the device (in Joules).

  • coulomb_peak_pos (1D array) – Coulomb peak positions for this device with respect to the reference gate potential (in volts).

  • zeeman (bool) – Whether or not Zeeman effects should be considered for the device

  • orb_B_effects (bool) – Whether or not orbital magnetic effects should be considered for the device

  • soc (3D array) – Describes the spin-orbit coupling in the system. The first index runs from 0 to 2 and describes coupling to the kx, ky, and kz operator respectively. The second and third indices run over the number of bands in the model. They define the matrices that are proportional to kx, ky, and kz for the specific type of spin-orbit coupling being considered.

Note

Default material values are for silicon. Electron affinity and bandgap are taken from http://www.virginiasemi.com/pdf/generalpropertiessi62002.pdf. Default donor (acceptor) binding energies are for phosphorus (boron).

__init__(mesh=None, inp_dict=None, conf_carriers='unspecified')

Constructor of the device class.

Parameters
  • mesh (mesh object, optional) – Mesh on which the device is defined. Default mesh if unspecified.

  • inp_dict (dict, optional) – Dictionary of device attributes. Every unspecified device attribute takes a default value instead.

  • conf_carriers (str or None, optional) – The type of charge carrier to be modeled quantum mechanically. Electrons (‘e’), holes (‘h’), or None’. Default: ‘unspecified’, in which case confined carriers are set to be electrons and a warning is raised.

add_to_V(V, region=None)

Add a constant potential to a given region

Parameters
  • V (float) – Potential to be added.

  • region (string or None) – If string, adds V to the potential defined over the region

add_to_ref_potential(phi, region=None)

Adds a constant potential to the reference potential in a given region.

Parameters
  • phi (float) – Potential to be added.

  • region (string or None, optional) – If string, adds phi to the reference potential defined over the region. If None, adds phi everywhere.

align_bands(ref_mat)

Sets the reference potential to be spatially dependent to have all bands aligned with respect to macroscopic potential parameters.

Parameters

ref_mat (material object) – Reference material.

Note

See Band alignment from atomistic simulations.

While the reference potential \(\phi_F\) is chosen to correctly align the bands, the reference potential \(\phi_0^\mathrm{ref}\) is set to have \(e\phi_0^\mathrm{ref}\) equal to the intrinsic work function of the reference material at zero temperature.

band_weight()

Computes the eigenfunction probability density to be in each band. Only works when there are multiple bands in the model (e.g. holes or electrons with spin).

Returns

2D array – Dimensions are number of eigenstates x bands. The entries contain the associated weights.

cond_band_edge()

Output conduction band edge at each mesh point.

Returns

2d Field – The conduction band edge in Joules at each local node.

get_N_particle_subspace(particle_number)

Gets the many-body subspace with a certain number of particles.

Parameters

particle_number (int) – The number of particles in the subspace to get.

Returns

Subspace object – A many-body subspace object.

get_Vconf()

Getter for the total confinement potential energy.

Returns

2d Field – The total confinement potential energy in Joules at each local node.

Note

The total confinement potential energy is Vconf = V + Vext.

get_Vmin()

Gets the minimum potential energy.

Returns

float – The minimum potential energy in Joules.

get_electron_g_factor()

Get the g-tensor for electrons

Returns

(2D or 3D array) – If the g-factor has been set by the user, a 2D array is returned. This array represents the isotropic g-tensor of the device. If the g-factor has not been set, it is constructed from the material parameters of the materials making up the device. In this case, a 3D array is returned, where the first index runs over the elements of the mesh (over which the device is defined) and the second and third indices run over directions.

Note

The g-tensor can be constructed from the g-factor because we assume that the tensor are isotropic (proportional to the identity matrix).

get_hole_Zeeman_params()

Get the g tensors for the holes.

Returns

(tuple) – The first entry of the tuple is the kappa g-tensor, and the second is the q g-tensor. If these have not been set by the user using the set_hole_zeeman_params<qtcad.device.device.Device.set_hole_zeeman_params() method, they are constructed from the materials defining the regions of the device and are position-dependent. The output tensors are 3D arrays in this case where the first index runs over the elements of the mesh and the second and third indices run over directions. In the case where the user sets the paramters values, using the set_hole_Zeeman_params method, they are assumed to be constant throughout the device and are returned as 2D arrays.

Note

The g-tensors can be constructed from g-factors because we assume

that tensors are isotropic (proportional to the identity matrix).

get_k_squared(state=0, directions=2)

Compute expectation value of momentum squared.

Parameters
  • state (int) – which state we wish to compute the expectation value of. Defaults to 0, i.e. the ground state

  • directions (list or int, optional) – list of directions for which we want the expectation values. 0 <-> ‘x’, i.e. <kx^2>, 1 <-> ‘y’, i.e. <ky^2>, and 2 <-> ‘z’, i.e. <kz^2>. If an integer is used, a single value will be output.

Returns

list or float – Expectation values of the required momentum operators squared.

Note

The device must contain a non-trivial eigenfunctions attribute in order for this method to be applied. This is typically achieved using the Schrodinger or Schrodinger-Poisson solver.

get_many_body_density(particle_number, eigenstate_index)

Gets the many-body particle density of a certain eigenstate for a given particle number.

Parameters
  • particle_number (int) – The number of particles in the many-body subspace of interest.

  • eigenstate_index (int) – The index of the many-body eigenstate within the subspace.

Returns

1D array – The particle density at each global node of the mesh.

get_potential_nrg()

Calculates the potential energy from the electric potential and band parameters.

load(path)

Loads a device from a json file.

new_dirichlet_bnd(bnd_spec, phi)

Create a Dirichlet boundary, i.e. an essential boundary in which the potential is set directly by the user (no built-in potential).

Parameters
  • bnd_spec (str or ndarray) – Boundary specifier.

  • phi (float) – Electric potential at the boundary.

Note

The boundary specifier is typically the string labeling the physical group for the boundary defined in the mesh file.

It may otherwise be a 1d array of integers labeling each global node in the boundary, or a 2d array of integers labeling each local node in the boundary.

new_frozen_bnd(bnd_spec, phi, material, doping, doping_type, binding_energy, dopant_degen=None)

Create a frozen boundary, i.e. an Ohmic boundary with a semiconductor at equilibrium, for which the electric potential is calculated assuming the chemical potential lies between the dopant level and the band edge, an assumption that only becomes valid at low temperature.

Parameters
  • bnd_spec (str or ndarray) – Boundary specifier.

  • phi (float) – Potential applied at the boundary.

  • material (material object) – The material for the semiconductor at the boundary.

  • doping (float) – The doping concentration of the semiconductor at the boundary.

  • doping_type (str) – The type of doping: “n” or “p”.

  • binding_energy (float) – The binding energy of the dopants

  • dopant_degen (int or None, optional) – The dopant level degeneracy. If None, uses default value: 2 for donors, 4 for acceptors.

Note

The boundary specifier is typically the string labeling the physical group for the boundary defined in the mesh file.

It may otherwise be a 1d array of integers labeling each global node in the boundary, or a 2d array of integers labeling each local node in the boundary.

new_gate_bnd(bnd_spec, phi, Ew)

Create a gate boundary, i.e. an essential boundary with a primary variable equal to the total potential due to the metal work function and an applied potential.

Parameters
  • bnd_spec (str or ndarray) – Boundary specifier.

  • phi (float) – Potential applied at the boundary.

  • Ew (float) – Metal work function.

Note

The boundary specifier is typically the string labeling the physical group for the boundary defined in the mesh file.

It may otherwise be a 1d array of integers labeling each global node in the boundary, or a 2d array of integers labeling each local node in the boundary.

new_insulator(bnd_nodes, set_no_charge=True)

Add an insulator boundary to the list of Schrodinger boundaries.

Parameters
  • bnd_nodes (1D array or str) – Indices of the global nodes that belong to the boundary, or region label (string).

  • set_no_charge (bool, optional) – If True, make sure that all classical charge vanishes inside the insulator nodes. Default: True.

new_ohmic_bnd(bnd_spec)

Create an Ohmic boundary, i.e. an essential boundary with a primary variable equal to the built-in potential determined by the the doping properties and material properties at the nodes belonging to the boundary.

Parameters

bnd_spec (str or ndarray) – Boundary specifier.

Note

The boundary specifier is typically the string labeling the physical group for the boundary defined in the mesh file.

It may otherwise be a 1d array of integers labeling each global node in the boundary, or a 2d array of integers labeling each local node in the boundary.

new_quantum_lead_bnd(bnd_spec, essential=False)

Adds a boundary interfacing with a quantum lead.

Parameters
  • bnd_spec (str or ndarray) – Boundary specifier.

  • essential (bool) – whether the boundary is essential (natural). essential=False if floating source (Neumann) boundary applies. In this case, the electrical boundary value is determined by resolving the electronic structure in the lead and its quantum statistics.

new_region(tag, material, pdoping=0, ndoping=0, ionization_model=None, Ea=None, Ed=None, gD=None, gA=None)

Create a region and add it to the device.

Parameters
  • tag (str or int) – The physical group tag (name or integer) to associate with the physical volume in the mesh file.

  • material (material object) – The material of which the region consists. This attribute is optional.

  • pdoping (float, optional) – The density of acceptors in cm^-3. Default: 0.

  • ndoping (float, optional) – The density of donors in cm^-3. Default: 0.

  • ionization_model (string, optional) – The ionization model used in this Region. Default: None, which means that the default ionization model stored in the device will be used, which itself is “complete” by default.

  • Ed (float, optional) – The donor binding energy. If None, does not modify the value currently stored in the device.

  • Ea (float, optional) – The acceptor binding energy. If None, does not modify the value currently stored in the device.

  • gD (int, optional) – The donor level degeneracy. If None, uses the default value for the current material.

  • gA (int, optional) – The acceptor level degeneracy. If None, uses the default value for the current material.

new_schottky_bnd(bnd_spec, phi, barrier)

Create a Schottky boundary, i.e. an essential boundary with a primary variable equal to the total potential due to the Schottky barrier height and an applied potential.

Parameters
  • bnd_spec (str or ndarray) – Boundary specifier.

  • phi (float) – Potential applied at the boundary.

  • barrier (float) – Schottky barrier height.

Note

The boundary specifier is typically the string labeling the physical group for the boundary defined in the mesh file.

It may otherwise be a 1d array of integers labeling each global node in the boundary, or a 2d array of integers labeling each local node in the boundary.

print_N_body_energies(N)

Prints the N-body energies.

Parameters

N (int) – Particle number of the subspace.

print_energies()

Print energy levels in eV.

Note

For holes the energies printed are those consistent with the convention defined in QTCAD convention for holes, i.e. the convention takes holes to be particles with positive mass and charge.

print_ith_evec(N, i)

Prints the energy and eigenvector of a given N-particles state.

Parameters
  • N (int) – Particle number of the subspace.

  • i (int) – eigenvalue index

reset_dot_region()

Resets this device in a state with no dot region defined.

save(path='device.hdf5', attrs=None)

Saves the device into a json or hdf5 file.

Parameters
  • path (string, optional) – The path to the file in which to save the device. Default: device.hdf5 in current dir.

  • attrs (list of str or None, optional) – Attributes to save in the json or hdf5 file. If None, save all of them.

  • bnd_list_names (list of str, optional) – Names of the device attributes that are actually lists of boundaries.

save_energies(path)

Save energies and population factors resulting from the self-consistent Schrodinger-Poisson solution in a text file.

set_Bfield(B, r0=(0, 0, 0), gauge='symmetric')

Setter for magnetic field.

Parameters
  • B (1D or 2D array) – Magnetic field defined over the device. If 1D array is used, the magnetic field is assumed to be constant throughout the device. In this case, a vector potential can be computed and orbital magnetic effects can be included in the model. Alternatively the magnetic field can be given at every global node of the mesh. In this case, we expect a 2D array where one index runs over the global nodes, the other runs over the directions of the magnetic field (indices are interchangeable).

  • r0 (tuple, optional) – Defaults to (0, 0, 0). Three-tuple defining the origin of the coordinate system. Used to shift the when computing the vector potential.

  • gauge (string, optional) – Defaults to “symmetric”. The gauge of the vector potential. Currently available: “symmetric” and “landau”

Note

The arguments r0 and gauge are only used when orbital magnetic

effects are modelled.

The orbital magnetic effects can only be modelled if the magnetic

field is constant.

Zeeman and/or orbital magnetic field effects can be turned off by

setting the zeeman and/or orb_B_effects attributes of to False

set_V(V, region=None)

Setter for the intrinsic confinement potential energy. Sets the this potential energy to arbitrary user-defined values without accounting for conduction and valence band edges.

Parameters
  • V (1D array, 2D array, function or float) – Value of the confinement potential at each node.

  • region (string or None) – If string, sets the value of the potential to V for every element in the region. In this case, V must be a float.

Note

The total confinement potential energy is Vconf = V + Vext.

set_V_from_phi()

Calculates the potential energy from the electric potential and band parameters.

set_Vext(Vext, region=None)

Setter for the external confinement potential energy

Parameters
  • Vext (1D array, 2D array, function or float) – Value of the external confinement potential at each node.

  • region (string or None) – If string, sets the value of the potential to Vext for every element in the region. In this case, Vext must be a float.

Note

The total confinement potential energy is Vconf = V + Vext.

set_Zeeman(value)

Setter for the zeeman attribute. Turns on/off Zeeman effect when solving the Schrodinger equation.

Parameters

value (boolean) – value of the zeeman attribute

set_acceptor_binding_nrg(Ea, region=None)

Sets the acceptor binding energy Ea.

Parameters
  • Ea (1D array, 2D array, function or float) – The values of the binding energy.

  • region (string or None, optional) – If string, sets the value of the binding energy to Ea for every element in this region. In this case, Ea must be a float.

set_acceptor_degen(gA, region=None)

Sets the acceptor level degeneracy gA.

Parameters
  • gA (1D array, 2D array, function or int) – The values of the binding energy.

  • region (string or None, optional) – If string, sets the value of the acceptor degeneracy to gA for every element in this region. In this case, gA must be a float.

set_acceptor_density(NA, region=None)

Sets the acceptor density NA.

Parameters
  • NA (1D array, 2D array, function or float) – The values of the acceptor density.

  • region (string or None, optional) – If string, sets the value of the acceptor density to NA for every element in this region. In this case, NA must be a float.

set_applied_potential(bnd_label, phi_a)

Set the applied potential at some boundary.

Parameters
  • bnd_label (str) – The boundary label

  • phi_a (float) – The potential to apply

Note

The boundary must support applied potentials. This is the case for, e.g., gate and Schottky boundaries.

set_charge_density(rho, region=None)

Sets the charge density. Used only by the linear Poisson solver.

Parameters
  • rho (1D array, 2D array, function or float) – Value of the charge density at each node.

  • region (string or None, optional) – If string, sets the value of the charge density to rho for every element in this region. In this case, rho must be a float.

Note

In the non-linear Poisson solver, the device attribute rho will be overwritten at every self-consistent iteration by the equilibrium semiconductor charges (electron and hole densities, and ionized donor and acceptor densities). To add a fixed background volume charge that will be accounted for by the non-linear Poisson solver, please see set_vol_charge_density.

set_cond_band_edge(EC)

Set the electric potential to correspond to a given conduction band edge.

Parameters

EC (1d array or 2d array) – The conduction band edge at each global or local node.

set_coulomb_mat(coulomb_mat)

Sets the Coulomb interaction matrix

Parameters

coulomb_mat (2d or 4d array) – The Coulomb interaction matrix. If a 2D array is used, it is assumed that the overlaps between basis states are neglected. Otherwise, the full 4D interaction matrix is used.

set_donor_binding_nrg(Ed, region=None)

Sets the donor binding energy Ed.

Parameters
  • Ed (1D array, 2D array, function or float) – The values of the binding energy.

  • region (string or None, optional) – If string, sets the value of the binding energy to Ed for every element in this region. In this case, Ed must be a float.

set_donor_degen(gD, region=None)

Sets the donor level degeneracy gD.

Parameters
  • gD (1D array, 2D array, function or int) – The values of the binding energy.

  • region (string or None, optional) – If string, sets the value of the donor degeneracy to gD for every element in this region. In this case, gD must be a float.

set_donor_density(ND, region=None)

Sets the donor density ND.

Parameters
  • ND (1D array, 2D array, function or float) – The values of the donor density.

  • region (string or None, optional) – If string, sets the value of the donor density to ND for every element in this region. In this case, ND must be a float.

set_dot_region(dot_region)

Define global nodes corresponding to the quantum dot (quantum confined) region.

Parameters

submesh (string, SubMesh or list) – Submesh or list of submeshes corresponding to the dot region, or string or list of strings labeling region of the device (Gmsh physical groups), or list of lists/tuples defining corners of a box for each dot

Note

The dot region is internally defined as a list of region labels or SubMesh objects. If this method is called multiple times on the same Device object, inputs are concatenated together in the list. For example, calling this method first with dot_region="region_1", and then with dot_region="region_2" results in a dot region which is the union of global nodes belonging to "region_1" and "region_2".

Classical charge density for confined carriers is set to zero inside the dot region when using the non-linear Poisson solver. Both classical and confined charge for electrons and holes (i.e., all charge except those due to ionized dopant densities) are set to zero inside the dot region when using the quantum-well solver.

set_fermi_lvl(EF, region=None)

Sets the Fermi energy level.

Parameters
  • EF (1D array, 2D array, function or float) – The values of the Fermi energy.

  • region (string or None, optional) – If string, sets the value of the Fermi energy to EF for every element in this region. In this case, EF must be a float.

set_g_star(g_star)

Setter for effective g-factor for electrons

Parameters

g_star (float) – g-factor associated with the coupling between the magnetic field and the angular momentum operators

Note

If this setter is used, the same isotropic g-tensor is assumed to

be valid for all nodes of the device.

set_hole_zeeman_params(kappa, q)

Setter for hole Zeeman parameters

Parameters
  • kappa (float) – g-factor associated with the coupling between the magnetic field and the angular momentum operators

  • q (float) – g-factor associated with the coupling between the magnetic field and the angular momentum operators to the third power

  • Note

    If this setter is used, the same isotropic g tensors are

    assumed to be valid for all nodes of the device.

set_ionization_model(model)

Sets the default ionization model of the device.

Parameters

model (string) – The model to be used.

Note

Available models are “complete” and “incomplete”.

set_orb_B_effects(value)

Setter for the orb_B_effects attribute. Turns on/off orbital magnetic effects when solving the Schrodinger equation.

Parameters

value (boolean) – value of the orb_B_effects attribute

set_permittivity(eps, region=None)

Sets the permittivity.

Parameters
  • eps (1D array, 2D array, function or float) – Value of the permittivity at each node.

  • region (string or None) – If string, sets the value of the permittivity to eps for every element in the region. In this case, eps must be a float.

set_potential(phi)

Set the electric potential.

Parameters

phi (1d array) – Value of the potential at each global node.

set_soc(soc)

Add a spin-orbit coupling term to the model.

Parameters

soc (list of 2D numpy arrays) – Each entry of the list contains the matrix proportional to kx, ky, and kz respectively which encodes a particular form of soc. If there is no coupling for a particular k direction, the entry in the tuple can be set to None

set_strain(strain, region=None)

Setter for strain throughout the device

Parameters
  • strain (numpy array or func which outputs an array) – 3 x 3 numpy array that defines the uniform strain throughout the device or a function of 3 variables (x, y, z) that outputs a 3 x 3 array defining the strain at every global node of the mesh.

  • region (string or None) – If string, sets the value of the attribute to val for every element in the region

set_surf_charge_density(surf_spec, density)

Sets the surface charge density at a boundary or interface.

Parameters
  • surf_spec (str or ndarray) – Surface specifier.

  • sigma (float) – Surface charge density.

Note

The surface specifier is typically the string labeling the physical group for the surface defined in the mesh file.

It may otherwise be a 1d array of integers labeling each global node in the surface, or a 2d array of integers labeling each local node in the surface.

set_temperature(T)

Sets the device temperature.

Parameters

T (float) – The device temperature, assumed uniform.

set_vlnce_band_edge(EV)

Set the electric potential to correspond to a given valence band edge.

Parameters

EV (1d array or 2d array) – The valence band edge at each global or local node.

set_vol_charge_density(rho_0, region=None)

Sets the background volume charge density.

Parameters
  • rho_0 (1D array, 2D array, function or float) – Value of the charge density at each node.

  • region (string or None, optional) – If string, sets the value of the background volume charge density to rho_0 for every element in this region. In this case, rho_0 must be a float.

Note

The background volume charge density attribute rho_0 is only used by the non-linear Poisson solver. To set the charge density in a way that will be accounted for by the linear Poisson solver, use set_charge_density instead, which directly sets the total charge density rho used by the linear Poisson solver.

update_mesh(mesh)

Update the mesh over which the device is defined.

Parameters

mesh (Mesh object) – The new mesh over which the device must now be defined.

vlnce_band_edge()

Output valence band edge at each mesh point.

Returns

2d Field – The conduction band edge in Joules at each local node.

warning_unassigned_region()

Issues a warning if there is any part of the device that is not assigned a region

class qtcad.device.device.SubDevice(parent, submesh)

Bases: device.core.device.SubDevice

A subdevice is a device built from a parent device from a submesh which is a subset of the mesh of the parent device.

__init__(parent, submesh)

Constructor of the SubDevice class.

Parameters
  • parent (Device object) – The Device object from which the subdevice is built.

  • submesh (SubMesh object) – The SubMesh object over which the subdevice is defined.

delete()

Removes the subdevice from memory.

Note

This also removes any reference to the subdevice from its parent device.