nanotools.relaxation.two_probe module

class nanotools.relaxation.two_probe.NBRelaxer(atomcell_L, atomcell_R, atomcell_C, transport_axis: str = 'Z', tolerance: float = 0.05, max_steps: int = 100, fix_2_layer_leads: bool = False, ase_atoms_L=None, ase_atoms_R=None, ase_atoms_C=None, atoms_C_relaxed=None)[source]

Bases: object

NBRelaxer class.

The NBRelaxer class is used for the relaxation of the center region of two-probe systems. It uses the pretrained universal neural network potential, CHGNet, to accurately predict atomic forces and energies, facilitating the structural relaxation process using the Atomic Simulation Environment (ASE).

atomcell_L

The AtomCell object representing the atomic structure of the left lead.

Type:

AtomCell

atomcell_R

The AtomCell object representing the atomic structure of the right lead.

Type:

AtomCell

atomcell_C

The AtomCell object representing the atomic structure of the central region.

Type:

AtomCell

transport_axis

The transport direction along which the relaxation is performed. Must be “X”, “Y”, or “Z”.

Type:

str

tolerance

The atomic force tolerance for the relaxation. The unit is eV/Å.

Type:

float

max_steps

The maximum number of optimization steps to perform.

Type:

int

fix_2_layer_leads

A boolean flag to fix more buffer atoms in the leads. Default is False.

Type:

bool

check_overlap()[source]

Check for overlap between the lead atoms and the buffer atoms in the center region.

Returns:

A tuple containing two lists:
  • list of indices of overlapping atoms between the leads and the buffer atoms in the center region

  • list of indices of non-overlapping atoms between the leads and the atoms in the center region

Return type:

tuple

Raises:

ValueError – If ase_atoms_C, ase_atoms_L, or ase_atoms_R is None.

check_overlap_side(atoms_C_side, atoms_side)[source]

Check for overlap between two sets of atoms.

Parameters:
  • atoms_C_side (Atoms) – First set of ASE Atoms.

  • atoms_side (Atoms) – Second set of ASE Atoms.

Returns:

Indices of atoms in atoms_side that overlap with any atom in atoms_C_side.

Return type:

ndarray

generate_pbc_for_direction(direction: int)[source]

Generate periodic boundary conditions based on the specified direction.

Parameters:

direction (int) – The direction along which the transport is happening. 0 for X, 1 for Y, 2 for Z.

Returns:

The periodic boundary conditions. A list of three booleans,

where the element at the index corresponding to the transport direction is False, and the others are True.

Return type:

list[bool]

Raises:

ValueError – If direction is not 0, 1, or 2.

parse_xyz_file(atom_cell, pbc=None)[source]

Parse the AtomCell object to an ASE Atoms object.

Parameters:
  • atom_cell (AtomCell) – An AtomCell object.

  • pbc (list of bool, optional) – The periodic boundary conditions. Defaults to [True, True, True].

Returns:

An ASE Atoms object, or None if atom_cell is None.

Return type:

Atoms

Raises:

ValueError – If pbc is not a list of three booleans.

plot_energy(filename: str | bytes | PathLike | None = None, show=False, trajectory_filename: str = 'relax.traj')[source]

Plot the energy during the optimization and save the plot to a file.

Parameters:
  • filename (os.Pathlike, optional) – The filename of the file to save the plot to. Defaults to “Energy”.

  • show (bool, optional) – Whether to display the plot. Defaults to False.

  • trajectory_filename (str) – The filename of the trajectory file to load.

Raises:

FileNotFoundError – If the trajectory file is not found.

Returns:

None

relax()[source]

Relax the structure.

This method converts AtomCell objects to ASE Atoms objects, sets the boundary condition of the center region to False along the transport direction, parses XYZ files, and relaxes the central region.

Raises:

ValueError – If atomcell_L, atomcell_R, or atomcell_C is None.

Returns:

None

relax_central_region()[source]

Relax the central region of the structure.

This method performs relaxation of the central region of the structure by setting up a calculator, fixing atoms, perturbating the structure, and performing relaxation using the BFGS minimizer of the ASE package.

Raises:

SystemExit – If there are insufficient buffer atoms in the center region.

Returns:

None

scale_atoms(atoms, scaling_factor, transport_direction)[source]

Scale the atoms along the specified transport direction.

Parameters:
  • atoms (Atoms) – ASE Atoms object.

  • scaling_factor (float) – The factor by which to scale the atoms.

  • transport_direction (int) – The direction along which to scale. Either 0 for X, 1 for Y, or 2 for Z.

Returns:

ASE Atoms object with scaled coordinates.

Return type:

Atoms

sort_atoms_along_direction(atoms: Atoms, direction_index: int)[source]

Sorts the atomic positions of an ASE Atoms object along a specified direction.

Parameters:
  • atoms (Atoms) – ASE Atoms object.

  • direction_index (int) – The direction along which to sort the atoms. Either 0, 1, or 2.

Returns:

ASE Atoms object with sorted atomic positions.

Return type:

Atoms

Raises:

ValueError – If direction_index is not 0, 1, or 2.

translate_atoms_along_direction(atoms: Atoms, direction: int, translation_distance: float)[source]

Translate the ASE Atoms object along the specified direction.

Parameters:
  • atoms (Atoms) – ASE Atoms object.

  • direction (int) – The direction along which to translate. Either 0 for X, 1 for Y, or 2 for Z.

  • translation_distance (float) – The distance to translate along the specified direction.

Returns:

ASE Atoms object with translated coordinates.

Return type:

Atoms

Raises:

ValueError – If direction is not 0, 1, or 2.

transport_direction_to_index(direction: str)[source]

Convert a transport direction (“X”, “Y”, or “Z”) to an index (0, 1, or 2).

Parameters:

direction (str) – The transport direction. Must be “X”, “Y”, or “Z”.

Returns:

The index corresponding to the transport direction. 0 for “X”, 1 for “Y”, and 2 for “Z”.

Return type:

int

Raises:

ValueError – If direction is not “X”, “Y”, or “Z”.

visualize_trajectory(trajectory_filename: str = 'relax.traj')[source]

Attempts to read a trajectory file and display it.

Parameters:

trajectory_filename (str) – The filename of the trajectory file to load.

Raises:

FileNotFoundError – If the trajectory file is not found, a warning message is printed and the program exits with status code 1.

Returns:

None

write(filename: str | bytes | PathLike | None = None, output_format: str = 'xyz')[source]

Write the relaxed atoms to a file in a specified format.

Parameters:
  • filename (os.Pathlike, optional) – The filename of the file. Defaults to “relaxed”.

  • output_format (str, optional) – The format to write the file in. Defaults to “xyz”.

Raises:

ValueError – If self.atoms_C_relaxed is None or if output_format is not a valid output format.

Returns:

None