qtcad.atoms.atoms module
Atomic structures.
- class Atoms(x: tuple | ndarray = None, y: tuple | ndarray = None, z: tuple | ndarray = None, mesh: Mesh = None, unit_cell: UnitCellZincblende = None, crystal_dir: ndarray = array([[1, 0, 0], [0, 1, 0], [0, 0, 1]]), pbc_directions: list = None, rng_seed: int = None, vff_params_name: str = 'Niquet', verbose: bool = True, legacy: bool = False)
Bases:
AtomsAtomic structure defined from layers and/or regions, with or without atomistic disorder in the form of random alloying and/or rough heterointerfaces, and with or without periodic boundary conditions.
- Attributes:
atom_num (int) – Number of atoms in the structure.
atom_species (Numpy.ndarray) – Vector whose elements are strings describing the atomic species (e.g. “Si”) of each atom in the structure.
atom_pos (Numpy.ndarray) – Matrix whose rows are the positions of the atoms in the structure; the ordering of the rows matches the ordering of atoms in
atom_species.atom_sublat (Numpy.ndarray) – Vector whose elements are integers describing the sublattice indices of each atom in the structure; the ordering of the elements matches the ordering of atoms in
atom_species. The sublattice indices may take values between 0 and 7, namely, there is one possible sublattice index for each atom in the conventional (cubic) unit cell.is_cation_site (Numpy.ndarray) – Vector whose elements are booleans specifying whether each atom in the structure is on a cation site (True) or on an anion site (False); the ordering of the elements matches the ordering of atoms in
atom_species. Used to define random alloys in III-V semiconductor structures; not physically meaningful in group-IV semiconductor structures.x (tuple or Numpy.ndarray) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the atomic structure along the \(x\) axis.
y (tuple or Numpy.ndarray) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the atomic structure along the \(y\) axis.
z (tuple or Numpy.ndarray) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the atomic structure along the \(z\) axis.
mesh (
Mesh) – 3D mesh whose geometry determines the shape of the atomic structure. Its physical volumes may define volumes within the atomic structure whose chemical composition may be altered.unit_cell (
UnitCellZincblende) – Unit cell to be sufficiently repeated to obtain an atomic structure filling the desired volume of space.crystal_dir (Numpy.ndarray) – Matrix of integers specifying crystal directions in terms of primitive lattice vectors of the conventional (cubic) unit cell. The first (second, third) row of the matrix specifies the crystal direction to be made parallel to the \(x\) (\(y\), \(z\)) axis.
pbc_directions (list of int) – List of integers chosen from 0, 1, and 2 specifying the directions (\(x\), \(y\), and \(z\), respectively) along which periodic boundary conditions are applied. If None, the atomic structure is finite along the \(x\), \(y\), and \(z\) axes. If not None, the atomic structure is periodic along the specified axis (or axes); the atomic structure is then the largest possible supercell that fits along this axis (or these axes) within the range (or ranges) specified through
x,y, and/orz.periods (Numpy.ndarray) – Vector containing the periods along the \(x\), \(y\) and \(z\) axes of the box. A period of infinity indicates that the atomic structure is non-periodic along the corresponding axis.
vff_params_name (str) – Name of the Keating valence force-field model parameter set, which is typically the name of the first author of the paper from which the parameters were extracted. Used to obtain the equilibrium bond lengths, which are used to establish which pairs of atoms are nearest neighbors as well as to compute the atomistic strain tensor, and in the Keating valence force-field model solver. The parameter set names include “Niquet” and “Gawarecki”. Set to “Niquet” unless specified through the
set_vff_paramsmethod.tb_params_name (str or tuple of str) – If a str, name of the tight-binding model parameter set, which is typically the name of the first author of the paper from which the parameters were extracted. If a tuple, names of the tight-binding model parameter sets ordered according to precedence; if multiple parameter sets are available for a given atom or pair of atoms, the parameter set is set to the available parameter set whose name first appears in the tuple. The parameter set names include “Niquet”, “Boykin”, “Vogl”, “Klimeck_cb”, “Klimeck_vb”, “Kim”, and “Gawarecki”. Only used in the atomistic tight-binding Schrödinger equation solver. Set to (“Niquet”, “Kim”) unless specified through the
set_tb_paramsmethod.verbose (bool) – Whether information pertaining to atomic structure manipulation is printed.
legacy (bool) – Whether to use the legacy algorithm for unit cell tiling, which is slower and significantly less memory-efficient than the new algorithm. The atomic structures resulting from the two algorithms may differ by a global shift of the positions of bulk atoms and in surface termination.
phi (callable) – Function which specifies the electric potential in the atomic structure. It takes as input three floats, the Cartesian coordinates of a point in space, and outputs a float, the value of the electric potential at that point. None (corresponding to zero electric potential) unless specified through the
set_potentialmethod. At most one ofphianddcan be specified at a time.d (
Device) – Device containing attributesmeshandphiwhich respectively describe the finite-element mesh and the electric potential obtained via a finite-element Poisson solver; used to specify the electric potential in the atomic structure. None (corresponding to zero electric potential) unless specified through theset_potential_from_devicemethod. At most one ofphianddcan be specified at a time.phi_to_add (list) – List of tuples specifying additional contributions to the electric potential in the atomic structure. Each tuple contains two elements: (1) a float or a function which specifies the additional contribution to the electric potential (floats define constant potentials and functions follow the same conventions as
phi) and (2) a string, an integer, or a list of strings and/or integers specifying the physical group tag(s) (name(s) or index/indices) of the physical volume(s) in the mesh file defining the atomic structure region over which to add the potential (or None to add the potential over the entire atomic structure). Empty list unless specified through theadd_to_potentialmethod.B (tuple) – Vector of the Cartesian coordinates of an external, constant magnetic field. None (corresponding to zero magnetic field) unless specified through the
set_Bfieldmethod.soc (bool) – Whether spin-orbit coupling is considered in the atomistic tight-binding Schrödinger equation solver. False unless specified through the
set_socmethod.zeeman (bool) – Whether the Zeeman effect is considered in the atomistic tight-binding Schrödinger equation solver. False unless specified through the
set_Zeemanmethod.energies (Numpy.ndarray) – Vector containing the eigenenergies obtained via the atomistic tight-binding Schrödinger equation solver, sorted in increasing order.
eigenfunctions (tuple of
Wavefunction) – The corresponding atomistic tight-binding eigenfunctions. The global phase of each eigenfunction is fixed so that the coefficient with the largest magnitude is real and positive.
- __init__(x: tuple | ndarray = None, y: tuple | ndarray = None, z: tuple | ndarray = None, mesh: Mesh = None, unit_cell: UnitCellZincblende = None, crystal_dir: ndarray = array([[1, 0, 0], [0, 1, 0], [0, 0, 1]]), pbc_directions: list = None, rng_seed: int = None, vff_params_name: str = 'Niquet', verbose: bool = True, legacy: bool = False) None
Check that the inputs are well-defined and tile the unit cell so as to fill the volume of atomic structure.
The geometry of the atomic structure may be specified by (1) the inputs
x,y, andz, in which case the atomic structure has the shape of a box (a rectangular parallelepiped) defined by these inputs, (2) the inputmesh, in which case the atomic structure has the shape of the mesh, or (3) the inputsx,y,z, andmesh, in which case the atomic structure has the shape of the intersection of the box and the mesh.- Parameters:
x (tuple or Numpy.ndarray, optional) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the atomic structure along the \(x\) axis. If None, it is determined from the bounding box of the input
mesh. Default: None.y (tuple or Numpy.ndarray, optional) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the atomic structure along the \(y\) axis. If None, it is determined from the bounding box of the input
mesh. Default: None.z (tuple or Numpy.ndarray, optional) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the atomic structure along the \(z\) axis. If None, it is determined from the bounding box of the input
mesh. Default: None.mesh (
Mesh, optional) – 3D mesh whose geometry determines the shape of the atomic structure. Its physical volumes may define volumes within the atomic structure whose chemical composition may be altered. Default: None.unit_cell (
UnitCellZincblende, optional) – Unit cell to be sufficiently repeated to obtain an atomic structure filling the desired volume of space. Default: unit cell of bulk, unstrained silicon.crystal_dir (Numpy.ndarray, optional) – Matrix of integers specifying crystal directions in terms of primitive lattice vectors of the conventional (cubic) unit cell. The first (second, third) row of the matrix specifies the crystal direction to be made parallel to the \(x\) (\(y\), \(z\)) axis. Default: identity matrix
Numpy.identity(3, dtype=int).pbc_directions (list of str, optional) – List of strings chosen from “x”, “y”, and “z” specifying the directions (\(x\), \(y\), and \(z\), respectively) along which periodic boundary conditions are applied. If None, the atomic structure is finite along the \(x\), \(y\), and \(z\) axes. If not None, the atomic structure is periodic along the specified axis (or axes); the atomic structure is then the largest possible supercell that fits along this axis (or these axes) within the range (or ranges) specified through
x,y, and/orz. Periodic boundary conditions are allowed only if the inputmeshis None. Default: None.rng_seed (int, optional) – Seed for the random number generator used to construct random alloys. If None, the seed is set from fresh, unpredictable entropy from the operating system. Default: None.
vff_params_name (str, optional) – Name of the Keating valence force-field model parameter set, which is typically the name of the first author of the paper from which the parameters were extracted. Used to obtain the equilibrium bond lengths, which are used to establish which pairs of atoms are nearest neighbors as well as to compute the atomistic strain tensor, and in the Keating valence force-field model solver. Surface atoms with fewer than two nearest neighbors are excluded from the atomic structure; the notion of nearest neighbor is established by the equilibrium bond lengths. The parameter set names include “Niquet” and “Gawarecki”. Default: “Niquet”.
verbose (bool, optional) – Whether information pertaining to atomic structure manipulation is printed. Default: True.
legacy (bool, optional) – Whether to use the legacy algorithm for unit cell tiling, which is slower and significantly less memory-efficient than the new algorithm. The atomic structures resulting from the two algorithms may differ by a global shift of the positions of bulk atoms and in surface termination. Default: False.
- add_to_potential(phi: float | Callable, tag: str | int | list | None = None) None
Add term to electric potential in a given region of the atomic structure.
This potential is used in the atomistic tight-binding Schrödinger equation solver, specifically by adding \(-e \varphi\) on the diagonal of the Hamiltonian, where \(e\) is the elementary charge and \(\varphi\) is the electric potential. The input potential is added to any potential specified via the
set_potentialorset_potential_from_devicemethods.- Parameters:
phi (float or callable) – Electric potential to add. If a float, the potential is assumed to be constant. If a function, it takes as input three floats, the Cartesian coordinates of a point in space, and outputs a float, the value of the electric potential at that point.
tag (str, int, or list, optional) – The physical group tag(s) (name, integer, or list of names and/or integers) of the physical volume(s) in the mesh file defining the region over which to add the potential. If None, the potential is added over the entire atomic structure. Default: None.
- get_atomistic_strain_tensor() ndarray
Compute atomistic strain tensor on each atomic site.
The implementation is based on Eq. 19 of [PKW+98]. In the reference unstrained atomic structure, the nearest neighbors of each atom are assumed to lie along the directions of the corners of a regular tetrahedron with bond lengths normalized to the equilibrium bond lengths.
- Returns:
Numpy.ndarray – Strain tensor for each atomic site; the first index is the atomic site index, and the second and third indices are Cartesian direction indices. The strain tensor is set to 0 for surface atoms.
- get_potential(method: str = 'linear') ndarray
Interpolate or evaluate the external electric potential on the positions of the atoms in the atomic structure.
- Parameters:
method (str, optional) – Method used to interpolate the finite-element electric potential on the atomic positions, if applicable. The options are “nearest” for a nearest-neighbor interpolation, “linear” for a linear interpolation, and “cubic” for a cubic interpolation. For interpolations from first-order finite-element meshes, “linear” offers optimal accuracy and speed. Default: “linear”.
- Returns:
Numpy.ndarray – Vector containing the values of the electric potential at the positions of the atoms.
- new_layer(atom_species: ndarray, method: str | None = 'stochastic', z_bot: float | RoughSurface | None = None, z_top: float | RoughSurface | None = None) None
Define the geometry and material composition of a heterostructure layer. Alter the material composition of the layer accordingly.
- Parameters:
atom_species (Numpy.ndarray) – For a heterostructure layer with a crystalline structure, vector whose elements are strings describing the atomic species (e.g. “Si”) of each atom of the conventional (cubic) unit cell in the layer; it should contain 8 elements whose group character (group III, group IV, group V) matches those of the elements of the attribute
unit_cell.atom_species_cc. For a randomly-alloyed heterostructure layer, matrix whose first column contains strings describing the atomic species (e.g. “Si”) of each atom in the random alloy, and whose second column contains floats or functions describing the ratio of the corresponding atomic species in the random alloy. Functions describing alloy ratios must take as input three floats (the Cartesian coordinates of a point in space) and output a float (the ratio of the corresponding atomic species at that point). Ratios must be between 0.0 and 1.0. For group-IV alloys, the ratios must sum to 1.0 at each point in space. For III-V alloys, the ratios of group-III (group-V) atomic species must sum to 1.0 at each point in space.method (str, optional) – For a randomly-alloyed heterostructure layer (i.e., for cases where
atom_specieshas two columns), method by which the random alloy is constructed. If “stochastic”, each atom in the heterostructure layer is independently replaced with one of the atomic species in the first column ofatom_speciesaccording to the probabilities provided in the second column; this may result in a chemical composition which slightly deviates from the target ratios, which reflects the inherent stochastic composition of finite random alloys. If “fixed”, the random alloy is generated such that the overall chemical composition exactly matches the ratios specified in the second column ofatom_species; this option is not available for ratios given as functions. This input is ignored ifatom_specieshas only one column. Default: “stochastic”.z_bot (float or
RoughSurface) – Bottom surface of the layer along the \(z\) axis. If a float, the surface is assumed to be perfectly flat. If None, the layer is assumed to extend to the minimal \(z\) position of the atomic structure. Default: None.z_top (float or
RoughSurface) – Top surface of the layer along the \(z\) axis. If a float, the surface is assumed to be perfectly flat. If None, the layer is assumed to extend to the maximal \(z\) position of the atomic structure. Default: None.
- new_region(atom_species: ndarray, method: str | None = 'stochastic', x_min: float | RoughSurface | None = None, x_max: float | RoughSurface | None = None, y_min: float | RoughSurface | None = None, y_max: float | RoughSurface | None = None, z_min: float | RoughSurface | None = None, z_max: float | RoughSurface | None = None, tag: str | int | list | None = None, is_inside: Callable | None = None) None
Define the geometry and material composition of a region. Alter the material composition of the region accordingly.
The geometry of the region may be specified by (1) the inputs
x_min,x_max,y_min,y_max,z_min, and/orz_max, in which case the region has a box-like shape (a rectangular parallelepiped whose faces may or may not be rough) defined by these inputs, (2) the inputtag, in which case the region is the volume with physical tagtagof the mesh file defining the atomic structure, (3) the inputis_inside, in which case the region is the volume of space for which the functionis_insidereturns True, or (4) any combination of the above inputs, in which case the region is the intersection of the volumes defined by each of the considered inputs.- Parameters:
atom_species (Numpy.ndarray) – For a region with a crystalline structure, vector whose elements are strings describing the atomic species (e.g. “Si”) of each atom of the conventional (cubic) unit cell in the region; it should contain 8 elements whose group character (group III, group IV, group V) matches those of the elements of the attribute
unit_cell.atom_species_cc. For a randomly-alloyed region, matrix whose first column contains strings describing the atomic species (e.g. “Si”) of each atom in the random alloy, and whose second column contains floats or functions describing the ratio of the corresponding atomic species in the random alloy. Functions describing alloy ratios must take as input three floats (the Cartesian coordinates of a point in space) and output a float (the ratio of the corresponding atomic species at that point). Ratios must be between 0.0 and 1.0. For group-IV alloys, the ratios must sum to 1.0 at each point in space. For III-V alloys, the ratios of group-III (group-V) atomic species must sum to 1.0 at each point in space.method (str, optional) – For a randomly-alloyed region (i.e., for cases where
atom_specieshas two columns), method by which the random alloy is constructed. If “stochastic”, each atom in the region is independently replaced with one of the atomic species in the first column ofatom_speciesaccording to the probabilities provided in the second column; this may result in a chemical composition which slightly deviates from the target ratios, which reflects the inherent stochastic composition of finite random alloys. If “fixed”, the random alloy is generated such that the overall chemical composition exactly matches the ratios specified in the second column ofatom_species; this option is not available for ratios given as functions. This input is ignored ifatom_specieshas only one column. Default: “stochastic”.x_min (float or
RoughSurface, optional) – Bottom surface of the region along the \(x\) axis. If a float, the surface is assumed to be perfectly flat. If None, the region is assumed to extend to the minimal \(x\) position of the atomic structure. Default: None.x_max (float or
RoughSurface, optional) – Top surface of the region along the \(x\) axis. If a float, the surface is assumed to be perfectly flat. If None, the region is assumed to extend to the maximal \(x\) position of the atomic structure. Default: None.y_min (float or
RoughSurface, optional) – Bottom surface of the region along the \(y\) axis. If a float, the surface is assumed to be perfectly flat. If None, the region is assumed to extend to the minimal \(y\) position of the atomic structure. Default: None.y_max (float or
RoughSurface, optional) – Top surface of the region along the \(y\) axis. If a float, the surface is assumed to be perfectly flat. If None, the region is assumed to extend to the maximal \(y\) position of the atomic structure. Default: None.z_min (float or
RoughSurface, optional) – Bottom surface of the region along the \(z\) axis. If a float, the surface is assumed to be perfectly flat. If None, the region is assumed to extend to the minimal \(z\) position of the atomic structure. Default: None.z_max (float or
RoughSurface, optional) – Top surface of the region along the \(z\) axis. If a float, the surface is assumed to be perfectly flat. If None, the region is assumed to extend to the maximal \(z\) position of the atomic structure. Default: None.tag (str, int, or list, optional) – The physical group tag(s) (name, integer, or list of names and/or integers) of the physical volume(s) in the mesh file defining the region. Default: None.
is_inside (callable, optional) – Function which specifies which atomic positions are inside the region. It takes as input three floats, the Cartesian coordinates of a point in space, and outputs a bool which is True if the point lies inside the region. Default: None.
- print_energies(decimals: int = 8) None
Print the eigenenergies obtained via the atomistic tight-binding Schrödinger equation solver in units of electronvolts.
- Parameters:
decimals (int, optional) – Number of decimals to print.
- save(path: str) None
Save the atomic structure to a file.
- Parameters:
path (str) – Path to the file in which to save the atomic structure; all file extensions are accepted.
- save_energies(path: str) None
Save to a
.txtfile the eigenenergies obtained via the atomistic tight-binding Schrödinger equation solver in units of electronvolts.- Parameters:
path (str) – Path to the
.txtfile in which to save the eigenenergies.
- save_xyz(path: str, comment: str = '') None
Save
.xyzfile of the atomic structure.- Parameters:
path (str) – Path in which the
.xyzfile is to be saved.comment (str, optional) – Content of the comment line of the
.xyzfile. Default:"".
Note
By convention, the atomic coordinates stored in the
.xyzfile are expressed in units of angstroms rather than metres.
- set_Bfield(B: tuple | ndarray) None
Set an external, constant magnetic field in the atomic structure. It is not possible to set a magnetic field in atomic structures with periodic boundary conditions.
This magnetic field is used in the atomistic tight-binding Schrödinger equation solver. The impact of the magnetic field on the offsite blocks of the tight-binding Hamiltonian is modeled through the Peierls substitution, whereby the standard momentum \(\mathbf{p}\) is effectively replaced by the generalized momentum \(\mathbf{p}-e\mathbf{A}\), where \(e\) is the elementary charge and \(\mathbf{A}\) is the magnetic vector potential. The impact of the magnetic field on the onsite blocks of the tight-binding Hamiltonian is (optionally) modeled through the Zeeman effect, whereby the Zeeman Hamiltonian \(-\mathbf{\mu}\cdot\mathbf{B}\), where \(\mathbf{\mu}\) is the magnetic moment and \(\mathbf{B}\) is the magnetic field, is effectively added to the Hamiltonian.
- Parameters:
B (tuple or Numpy.ndarray) – Vector of the Cartesian coordinates of an external, constant magnetic field.
- set_Zeeman(value: bool) None
Set whether the Zeeman effect is considered in the atomistic tight-binding Schrödinger equation solver.
- Parameters:
value (bool) – Whether the Zeeman effect is considered.
- set_potential(phi: Callable) None
Set the electric potential in the atomic structure from a function.
This potential is used in the atomistic tight-binding Schrödinger equation solver, specifically by adding \(-e \varphi\) on the diagonal of the Hamiltonian, where \(e\) is the elementary charge and \(\varphi\) is the electric potential.
- Parameters:
phi (callable) – Function which specifies the electric potential in the atomic structure. It takes as input three floats, the Cartesian coordinates of a point in space, and outputs a float, the value of the electric potential at that point.
Note
This method sets the attribute
dto None.
- set_potential_from_device(d: Device) None
Set the electric potential in the atomic structure from the results of a finite-element simulation, taking into account the periodic boundary conditions that may exist in the device.
This potential is used in the atomistic tight-binding Schrödinger equation solver, specifically by adding \(-e (\varphi - \Delta\varphi_F)\) on the diagonal of the Hamiltonian, where \(e\) is the elementary charge, \(\varphi\) is the electric potential, and \(\Delta\varphi_F\) accounts for the difference in the gauges assumed in finite-element and atomistic simulations. Specifically, \(e \Delta\varphi_F = E_w^{\mathrm{Si}}(T) - E_w^{\mathrm{Si}}(0)\), where \(E_w^{\mathrm{Si}}(T)\) (\(E_w^{\mathrm{Si}}(0)\)) is the workfunction of intrinsic silicon at the device temperature \(T\) (at absolute zero).
- Parameters:
d (
Device) – Device containing attributesmeshandphiwhich respectively describe the finite-element mesh and the electric potential obtained via a finite-element Poisson solver, and on which thenew_periodic_bndmethod has been called to define periodic boundary conditions, if applicable.
Note
This method sets the attribute
phito None.
- set_soc(value: bool) None
Set whether spin-orbit coupling is considered in the atomistic tight-binding Schrödinger equation solver.
- Parameters:
value (bool) – Whether spin-orbit coupling is considered.
- set_tb_params(tb_params_name: str | Tuple[str, ...]) None
Set the parameter set for the atomistic tight-binding Schrödinger equation solver.
- Parameters:
tb_params_name (str or tuple of str) – If a str, name of the tight-binding model parameter set, which is typically the name of the first author of the paper from which the parameters were extracted. If a tuple, names of the tight-binding model parameter sets ordered according to precedence; if multiple parameter sets are available for a given atom or pair of atoms, the parameter set is set to the available parameter set whose name first appears in the tuple. The parameter set names include “Niquet”, “Boykin”, “Vogl”, “Klimeck_cb”, “Klimeck_vb”, “Kim”, and “Gawarecki”.
Note
The parameter sets are defined in the
materialsmodule.
- set_vff_params(vff_params_name: str) None
Set the parameter set to compute the atomistic strain tensor and for the Keating valence force-field model solver.
- Parameters:
vff_params_name (str) – Name of the Keating valence force-field model parameter set, which is typically the name of the first author of the paper from which the parameters were extracted. The parameter set names include “Niquet” and “Gawarecki”.
Note
The parameter sets are defined in the
materialsmodule.
- strain_atomic_structure(strain_tensor: Callable) None
Apply strain to atomic structure. This method is only available for atomic structures without periodic boundary conditions.
The displacement of the atoms from their current position is obtained by integrating the strain tensor along a line joining an arbitrary reference point to the position of the atom. The resulting atomic positions are stored in
atom_pos.- Parameters:
strain_tensor (callable) – Function which specifies the desired strain tensor in the atomic structure. It takes as input three floats, the Cartesian coordinates of a point in space, and outputs a symmetric \(3\times 3\) matrix (Numpy.ndarray), the value of the strain tensor at that point.
Note
Since this method is based on the continuum, infinitesimal strain tensor formalism, it may result in atomic structures whose atomistic strain tensor (as computed from
get_atomistic_strain_tensor) deviates from the input strain tensor.
- class SubAtoms(parent: Atoms, x: tuple | ndarray = None, y: tuple | ndarray = None, z: tuple | ndarray = None, tag: str | int | list = None, is_inside: Callable = None, pbc_directions: list = None, remove_atoms_with_few_nn: bool | None = True)
Bases:
SubAtomsRestricted atomic structure (typically used as an input to the atomistic tight-binding Schrödinger equation solver).
- __init__(parent: Atoms, x: tuple | ndarray = None, y: tuple | ndarray = None, z: tuple | ndarray = None, tag: str | int | list = None, is_inside: Callable = None, pbc_directions: list = None, remove_atoms_with_few_nn: bool | None = True) None
Check that the inputs are well-defined and restrict the parent atomic structure to the specified volume of space.
The geometry of the restricted atomic structure may be specified by (1) the inputs
x,y, andz, in which case the restricted atomic structure is the intersection of the parent atomic structure and a box (a rectangular parallelepiped) defined by these inputs, (2) the inputtag, in which case the restricted atomic structure is the intersection of the parent atomic structure and the volume with physical tagtagof the mesh of the parent atomic structure, (3) the inputis_inside, in which case the restricted atomic structure is the intersection of the parent atomic structure and the volume of space for which the functionis_insidereturns True, or (4) any combination of the above inputs, in which case the restricted atomic structure is the intersection of the parent atomic structure and the volumes defined by each of the considered inputs.- Parameters:
parent (
Atoms) – Atomic structure from which to extract the restricted atomic structure.x (tuple or Numpy.ndarray, optional) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the restricted atomic structure along the \(x\) axis. If None, it is automatically determined from the other inputs. This range may not be narrower than that of the parent atomic structure if periodic boundary conditions are applied along the \(x\) axis. Default: None.
y (tuple or Numpy.ndarray, optional) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the restricted atomic structure along the \(y\) axis. If None, it is automatically determined from the other inputs. This range may not be narrower than that of the parent atomic structure if periodic boundary conditions are applied along the \(y\) axis. Default: None.
z (tuple or Numpy.ndarray, optional) – Vector with two elements specifying the minimal and maximal positions of the bounding box of the restricted atomic structure along the \(z\) axis. If None, it is automatically determined from the other inputs. This range may not be narrower than that of the parent atomic structure if periodic boundary conditions are applied along the \(z\) axis. Default: None.
tag (str, int, or list, optional) – The physical group tag(s) (name, integer, or list of names and/or integers) of the physical volume(s) in the mesh file defining the volume of the restricted atomic structure. Default: None.
is_inside (callable, optional) – Function which specifies which atomic positions are inside the restricted atomic structure. It takes as input three floats, the Cartesian coordinates of a point in space, and outputs a bool which is True if the point lies inside the restricted atomic structure. Must be None if periodic boundary conditions are considered. Default: None.
pbc_directions (list of str, optional) – List of strings chosen from “x”, “y”, and “z” specifying the directions (\(x\), \(y\), and \(z\), respectively) along which periodic boundary conditions are applied. If the empty list, the atomic structure is finite along the \(x\), \(y\), and \(z\) axes. If None, the periodicity is inherited from the parent atomic structure. Periodic boundary conditions along a given direction are allowed if and only if the parent atomic structure has periodic boundary conditions along this direction. Default: None.
remove_atoms_with_few_nn (bool, optional) – Whether to remove surface atoms that have fewer than two nearest neighbors. Should only be set to False for visualization purposes. Default: True.