tbplas.Sample

class tbplas.Sample(*args: Union[SuperCell, SCInterHopping])

Interface class to FORTRAN backend.

_sc_list

list of supercells within the sample

Type:

List[SuperCell]

_hop_list

list of inter-hopping sets between supercells within the sample

Type:

List[SCInterHopping]

orb_eng

on-site energies of orbitals in the supercell in eV

Type:

(num_orb_tot,) float64 array

orb_pos

Cartesian coordinates of orbitals in the supercell in nm

Type:

(num_orb_tot, 3) float64 array

hop_i

row indices of hopping terms reduced by conjugate relation

Type:

(num_hop_tot,) int64 array

hop_j

column indices of hopping terms reduced by conjugate relation

Type:

(num_hop_tot,) int64 array

hop_v

energies of hopping terms in accordance with hop_i and hop_j in eV

Type:

(num_hop_tot,) complex128 array

dr

distances of hopping terms in accordance with hop_i and hop_j in nm

Type:

(num_hop_tot, 3) float64 array

_rescale

rescaling factor for the Hamiltonian reserved for compatibility with old version of TBPlaS

Type:

float

Notes

If periodic conditions are enabled, orbital indices in hop_j may be wrapped back if it falls out of the supercell. Nevertheless, the distances in dr are still the ones before wrapping. This is essential for adding magnetic field, calculating band structure and many properties involving dx and dy.

__init__(*args: Union[SuperCell, SCInterHopping]) None
Parameters:

args – supercells and inter-hopping sets within this sample

Returns:

None

Raises:
  • SampleVoidError – if len(args) == 0

  • SampleCompError – if any argument in args is not instance of ‘SuperCell’ or ‘SCInterHopping’ classes

  • SampleClosureError – if any ‘SCInterHopping’ instance has supercells not included in the sample

Methods

__init__(*args)

param args:

supercells and inter-hopping sets within this sample

add_subscriber(sub_name, sub_obj)

Add a new subscriber.

build_dxy_csr()

Build sparse dx and dy matrices in CSR format for TESTING purposes.

build_ham_csr()

Build sparse Hamiltonian for DOS and LDOS calculations using TBPM.

build_ham_dxy([orb_eng_cutoff, algo, sort_col])

Build the arrays for conductivity calculations using TBPM.

calc_bands(k_points[, enable_mpi, echo_details])

Calculate band structure along given k_path.

calc_dos(k_points[, enable_mpi, echo_details])

Calculate density of states for given energy range and step.

check_lock()

Check the lock state of the object.

init_hop([force_init])

Initialize hop_i, hop_j, hop_v, dr and reset rescale on demand.

init_orb_eng([force_init])

Initialize self.orb_eng on demand.

init_orb_pos([force_init])

Initialize self.orb_pos on demand.

load_array([data_dir])

Load array attributes from disk.

lock(sub_name)

Lock the object.

plot([fig_name, fig_size, fig_dpi, ...])

Plot lattice vectors, orbitals, and hopping terms.

rescale_ham([factor])

Rescale orbital energies and hopping terms.

reset_array()

Reset all modifications to self._orb_*, self._hop_* and self.dr.

save_array([data_dir])

Save array attributes and scaling factor to disk.

set_ham_csr(k_point[, convention])

Set up sparse Hamiltonian in csr format for given k-point.

set_ham_dense(k_point, ham_dense[, convention])

Set up Hamiltonian for given k-point.

set_k_point(k_point)

Set the k-point of Hamiltonian.

set_magnetic_field(intensity[, gauge])

Apply magnetic field perpendicular to xOy-plane to -z direction via Peierls substitution.

sync_array(**kwargs)

Abstract interface for 'sync_array' method of derived classes.

unlock()

Unlock the object.

update()

Interface for keeping data consistency in observer pattern.

Attributes

area_unit_cell

Get the area formed by a1 and a2 of the primitive cell.

energy_range

Get energy range to consider in calculations.

extended

Get the number of extended times of primitive cell.

nr_orbitals

Get the number of orbitals of primitive cell.

num_hop

Get the total number of hopping terms of the sample.

num_orb

Get the total number of orbitals of the sample.

rescale

Interface the '_rescale' attribute.

sc0

Interface for the 1st supercell.

volume_unit_cell

Get the volume of primitive cell.

__init__(*args: Union[SuperCell, SCInterHopping]) None
Parameters:

args – supercells and inter-hopping sets within this sample

Returns:

None

Raises:
  • SampleVoidError – if len(args) == 0

  • SampleCompError – if any argument in args is not instance of ‘SuperCell’ or ‘SCInterHopping’ classes

  • SampleClosureError – if any ‘SCInterHopping’ instance has supercells not included in the sample

add_subscriber(sub_name: str, sub_obj: Any) None

Add a new subscriber.

Parameters:
  • sub_name – subscriber name

  • sub_obj – subscriber object

Returns:

None

build_dxy_csr() Tuple[csr_matrix, csr_matrix]

Build sparse dx and dy matrices in CSR format for TESTING purposes.

NOTE: As zero elements are removed automatically when creating sparse matrices using scipy, Hamiltonian and dx/dy may have different indices and indptr. So, DO NOT use this method to generate dx and dy for TBPM calculations. Use the ‘build_ham_dxy’ method instead.

Returns:

(dx_csr, dy_csr) sparse dx and dy matrices in CSR format

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

build_ham_csr() csr_matrix

Build sparse Hamiltonian for DOS and LDOS calculations using TBPM.

Returns:

sparse Hamiltonian matrix in CSR format

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

build_ham_dxy(orb_eng_cutoff: float = 1e-05, algo: str = 'fast', sort_col: bool = False) Tuple[ndarray, ndarray, ndarray, ndarray, ndarray]

Build the arrays for conductivity calculations using TBPM.

NOTE: Two algorithms are implemented to build the arrays. “fast” is more efficient and uses less memory, but relies heavily on delicate pointer operations. “safe” is built on numpy functions and considered to be more robust, but uses more memory and twice slower than “fast”. In short, use “fast” for production runs and “safe” for testing purposes.

Parameters:
  • orb_eng_cutoff – cutoff for orbital energies in eV Orbital energies below this value will be dropped when constructing sparse Hamiltonian.

  • algo – specifies which core function to call to generate the arrays Use “fast” for production runs and “safe” for testing purposes.

  • sort_col – whether to sort column indices and data of CSR matrices after creation, for TESTING purposes Sorting the columns will take much time, yet does not boost TBPM calculations. DO NOT enable this option for production runs.

Returns:

(indptr, indices, hop, dx, dy) indptr: int64 array indices: int64 array hop: complex128 array dx: float64 array dy: float64 array

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

calc_bands(k_points: ndarray, enable_mpi: bool = False, echo_details: bool = True, **kwargs) Tuple[ndarray, ndarray]

Calculate band structure along given k_path.

Parameters:
  • k_points – (num_kpt, 3) float64 array FRACTIONAL coordinates of the k-points

  • enable_mpi – whether to enable parallelization over k-points using mpi

  • echo_details – whether to output parallelization details

  • kwargs – arguments for ‘calc_bands’ of ‘DiagSolver’ class

Returns:

(k_len, bands) k_len: (num_kpt,) float64 array in 1/NM length of k-path in reciprocal space, for plotting band structure bands: (num_kpt, num_orb) float64 array Energies corresponding to k-points in eV

calc_dos(k_points: ndarray, enable_mpi: bool = False, echo_details: bool = True, **kwargs) Tuple[ndarray, ndarray]

Calculate density of states for given energy range and step.

Parameters:
  • k_points – (num_kpt, 3) float64 array FRACTIONAL coordinates of the k-points

  • enable_mpi – whether to enable parallelization over k-points using mpi

  • echo_details – whether to output parallelization details

  • kwargs – arguments for ‘calc_dos’ of ‘DiagSolver’ class

Returns:

(energies, dos) energies: (num_grid,) float64 array energy grid corresponding to e_min, e_max and e_step dos: (num_grid,) float64 array density of states in states/eV

check_lock() None

Check the lock state of the object.

Returns:

None

Raises:

LockError – if the object is locked

init_hop(force_init: bool = False) None

Initialize hop_i, hop_j, hop_v, dr and reset rescale on demand.

If the arrays are None, build them from scratch. Otherwise, build them only when force_init is True.

Parameters:

force_init – whether to force initializing the arrays from scratch even if they have already been initialized

Returns:

None

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

init_orb_eng(force_init: bool = False) None

Initialize self.orb_eng on demand.

If self.orb_eng is None, build it from scratch. Otherwise, build it only when force_init is True.

Parameters:

force_init – whether to force initializing the array from scratch even if it has already been initialized

Returns:

None

init_orb_pos(force_init: bool = False) None

Initialize self.orb_pos on demand.

If self.orb_pos is None, build it from scratch. Otherwise, build it only when force_init is True.

Parameters:

force_init – whether to force initializing the array from scratch even if it has already been initialized

Returns:

None

load_array(data_dir: str = 'sample') None

Load array attributes from disk.

Parameters:

data_dir – directory in which data are saved

Returns:

None

lock(sub_name: str) None

Lock the object. Modifications are not allowed then unless the ‘unlock’ method is called.

Parameters:

sub_name – identifier of the subscriber

Returns:

None

Raises:

ValueError – if sub_name is not in subscribers

plot(fig_name: str = None, fig_size: Tuple[float, float] = None, fig_dpi: int = 300, sc_orb_colors: List[Callable[[ndarray], List[str]]] = None, sc_hop_colors: List[str] = None, inter_hop_colors: List[str] = None, **kwargs) None

Plot lattice vectors, orbitals, and hopping terms.

If figure name is give, save the figure to file. Otherwise, show it on the screen.

Parameters:
  • fig_name – file name to which the figure will be saved

  • fig_size – size of the figure

  • fig_dpi – resolution of the figure file

  • sc_orb_colors

  • sc_hop_colors – colors for the hopping terms of each supercell

  • inter_hop_colors – colors for the hopping terms each inter-hopping container

  • kwargs – arguments for the ‘plot’ method of ‘Super’ and ‘SCInterHopping’ classes

Returns:

None

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

  • ValueError – if view is illegal

rescale_ham(factor: float = None) None

Rescale orbital energies and hopping terms.

Reserved for compatibility with old version of TBPlaS.

All orbital energies and hopping terms will be divided by the factor w.r.t. their original values in primitive cell. So only the last call to this method will take effect.

Choose the factor such that the absolute value of the largest eigenvalue is smaller than 1 after rescaling. If no value is chosen, a reasonable value will be estimated from the Hamiltonian.

Parameters:

factor – rescaling factor

Returns:

None

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

reset_array() None

Reset all modifications to self._orb_*, self._hop_* and self.dr.

Returns:

None

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

save_array(data_dir: str = 'sample') None

Save array attributes and scaling factor to disk.

Parameters:

data_dir – directory to which data will be saved

Returns:

None

set_ham_csr(k_point: ndarray, convention: int = 1) csr_matrix

Set up sparse Hamiltonian in csr format for given k-point.

This is the interface to be called by external exact solvers. The callers are responsible to call the ‘init_*’ method.

Parameters:
  • k_point – (3,) float64 array FRACTIONAL coordinate of the k-point

  • convention – convention for setting up the Hamiltonian

Returns:

sparse Hamiltonian

Raises:

ValueError – if convention != 1

set_ham_dense(k_point: ndarray, ham_dense: ndarray, convention: int = 1) None

Set up Hamiltonian for given k-point.

This is the interface to be called by external exact solvers. The callers are responsible to call the ‘init_*’ methods and initialize ham_dense as a zero matrix.

Parameters:
  • k_point – (3,) float64 array Fractional coordinate of the k-point

  • ham_dense – (num_orb, num_orb) complex128 array incoming Hamiltonian

  • convention – convention for setting up the Hamiltonian

Returns:

None

Raises:

ValueError – if convention != 1

set_k_point(k_point: ndarray) None

Set the k-point of Hamiltonian.

Parameters:

k_point – (3,) float64 array FRACTIONAL coordinate of the k-point

Returns:

None

set_magnetic_field(intensity: float, gauge: int = 0) None

Apply magnetic field perpendicular to xOy-plane to -z direction via Peierls substitution.

The gauge-invariant Peierls phase from R_i to R_j is evaluated as

exp[-1j * |e|/(2*h_bar*c) * (R_i - R_j) * (A_i + A_j)]

or

exp[1j * |e|/(2*h_bar*c) * (R_j - R_i) * (A_j + A_i)]

For the reference, see eqn. 10-11 of ref [1] and eqn. 19-20 of ref [2]. Be aware that ref [1] uses |e| while ref [2] uses e with sign.

For perpendicular magnetic pointing to +z, we have A = (-By, 0, 0). Then

(R_j - R_i) * (A_j + A_i) = -B * (y_j + y_i) * (x_j - x_i)

However, there isn’t the minus sign in the source code of TBPLaS. This is because we are considering a magnetic field pointing to -z.

Note that the above formulae use Gaussian units. For SI units, the speed of light vanishes. This is verified by checking the dimension under SI units:

[e]/[h_bar] = IT/(L^2 M T^-1) = L^-2 M^-1 I T^-2 [R][A] = L * (M T^-2 I^-1 L) = L^2 M I^-1 T^-2

which cancel out upon multiplication. Similarly, under Gaussian units we have: [e]/[c*h_bar] = L^(-3/2) M^(-1/2) T [R][A] = L^(3/2) M^(1/2) T which also cancel out upon multiplication.

The analysis also inspires to calculate the scaling factor in SI units as:

|e/2h_bar| = 759633724404755.2

However, we use nm for lengths. So the scaling factor should be divided by 1e18, yielding 0.0007596337244047553. That is exactly the magic number pi / 4135.666734 = 0.0007596339008078771 in the source code of TBPLaS.

Reference: [1] https://journals.aps.org/prb/pdf/10.1103/PhysRevB.51.4940 [2] https://journals.jps.jp/doi/full/10.7566/JPSJ.85.074709

Parameters:
  • intensity – magnetic B field in Tesla

  • gauge – gauge of vector potential of the magnetic field to -z 0 for (By, 0, 0), 1 for (0, -Bx, 0), 2 for (0.5By, -0.5Bx, 0)

Returns:

None

Raises:
  • InterHopVoidError – if any inter-hopping set is empty

  • ValueError – if duplicate terms have been detected in any inter-hopping

  • ValueError – if gauge is not in (0, 1, 2)

sync_array(**kwargs) None

Abstract interface for ‘sync_array’ method of derived classes.

unlock() None

Unlock the object. Modifications are allowed then.

Returns:

None

update() None

Interface for keeping data consistency in observer pattern.

Returns:

None

__weakref__

list of weak references to the object (if defined)

property area_unit_cell: float

Get the area formed by a1 and a2 of the primitive cell.

Reserved for compatibility with old version of TBPlaS.

Returns:

area formed by a1 and a2 in NM^2

property energy_range: float

Get energy range to consider in calculations.

Reserved for compatibility with old version of TBPlaS.

All eigenvalues are between (-en_range / 2, en_range / 2) in eV.

Returns:

the energy range

property extended: float

Get the number of extended times of primitive cell.

Reserved for compatibility with old version of TBPlaS.

Returns:

number of extended times of primitive cells

property nr_orbitals: int

Get the number of orbitals of primitive cell.

Reserved for compatibility with old version of TBPlaS.

Returns:

nr_orbitals: integer number of orbitals in the primitive cell.

property num_hop: int

Get the total number of hopping terms of the sample.

Returns:

total number of hopping terms

property num_orb: int

Get the total number of orbitals of the sample.

Returns:

total number of orbitals

property rescale: float

Interface the ‘_rescale’ attribute.

Returns:

the rescaling factor

property sc0: SuperCell

Interface for the 1st supercell.

Returns:

the 1st supercell

property volume_unit_cell: float

Get the volume of primitive cell.

Reserved for compatibility with old version of TBPlaS.

Returns:

volume of primitive cell in NM^3