Atomic Simulation Environment interface

Introduction

The viperleed.calc package of supports running directly from surface structures from the Atomic Simulation Environment (ASE). ASE is a widely used Python library for setting up and controlling atomic simulations such as density functional theory calculations and molecular dynamics simulations. Via the ASE interface provided in viperleed.calc, an ase.Atoms object can be used as surface model for a LEED-I(V) calculation.

This can be used to enable automated runs and sample surface structures generated by any external software. ViPErLEED can act as a “back-end” for high-throughput calculations which accepts (almost) arbitrary structure model and returns raw I(V) curves or the agreement (in terms of \(R\) factor) with an experimental set of beams.

Surface structures passed to viperleed.calc via the API need to fulfill the same conventions (e.g. vacuum side of surface towards +z direction) as applicable for POSCAR files.

The ASE API provides a number of Python function that allow calling and starting LEED-I(V) calculations. For an example of how to use the ASE interface, see this tutorial. Please note that the ASE interface may require installing the R-factor extension as described in the installation instructions.

Note

The ViPErLEED ASE interface is an area of active development and may change in the future. If you have any suggestions for features or feedback concerning the ASE interface, feel free to contact the ViPErLEED team (preferably by opening an issue on the ViPErLEED GitHub page).

Automatic site definitions

To facilitate batch processing of structures, the ASE interface allows for a simplistic automated assignment of site definitions (usually done via the SITE_DEF parameter).

If the SITE_DEF parameter is not defined in a ViPErLEED run using the ASE interface, ViPErLEED will try to assign “surface” sites on its own. For this calculation, every atom is considered as a solid sphere with a radius proportional to the elements’ covalent radius. Then, going from highest to lowest atom (\(z\) position), ViPErLEED checks if an atom is “visible” from vacuum, or if the line of sight is blocked by higher-up atoms. Atoms that are “visible” are declared as surface atoms (site surf), while all other atoms will be given the default site type (def). This will create, at maximum, two site types per species. If you want to define additional site types, it is currently necessary to define them explicitly in the POSCAR file.

A collection of functions that run ViPErLEED from a provided ASE object.

class viperleed.calc.from_ase.SlabTransform(orthogonal_matrix: Any | None = None, uc_scaling: Any | None = None, cut_cell_c_fraction: float = 0.0)[source]

A container for unit-cell transformations.

orthogonal_matrix

Shape (3, 3). Can contain an arbitrary orthogonal transform (i.e., a combination of mirroring and rotations) that will be applied to the cell. The transformation is applied BEFORE scaling, i.e., using the original unit cell. The matrix (O) is applied to both the unit cell and all Cartesian atomic positions. The unit cell U (unit vectors as columns) is transformed to O @ U. Atom Cartesian coordinates (v) are transformed to O @ v. This attribute can be used, e.g., to apply a rotation or a mirroring of the slab.

Type:

numpy.ndarray, optional

uc_scaling

Scaling to be applied to the cell. Scaling is applied AFTER the transformation, i.e., along the new unit cell vectors. This can be used to apply an isotropic scaling in all directions or to stretch/compress along unit cell vectors in order to change lattice constants in some direction. If the scaling is given as a number or a 1-element sequence, an isotropic scaling will be applied to the unit cell and atom positions. If a sequence with three entries is provided, the scaling will be applied along the unit cell vectors in the given order. Specifically, given uc_scaling == (s1, s2, s3), the new unit cell vectors will be a’ = a * s1, b’ = b * s2, c’ = c * s3. Default is None, i.e., no scaling.

Type:

float or Sequence of float, optional

cut_cell_c_fraction

Fractional position along the c vector to cut the slab. (ViPErLEED expects the cell to have the surface at the “top” – i.e., at high z coordinates, and bulk-like layers at the bottom). Cutting will be performed along the vector c AFTER all transformations (see uc_transformation_matrix) are applied. Only the part of the cell above (i.e., >=) this value will be kept. Everything else is discarded. Default is zero.

Type:

float, optional

viperleed.calc.from_ase.plot_iv_from_csv(beam_file, output_file=None, beam_file_is_content=False, which_beams=[], legends=None)[source]

Plots LEED-I(V) spectra directly from a CSV file.

This function can be used to plot LEED-I(V) spectra directly from a CSV file without running any other calculations. It is assumed that the CSV file uses the format used by ViPErLEED. Returns nothing but writes the pdf file.

Parameters:
  • beam_file (string) – If beam_file_is_content is set as True, the string is expected to contain the contents of the CSV file, otherwise it should be a path to the file.

  • beam_file_is_content (bool, default = False) – If set to True, it is assumed that the beam_file contains a string with the CSV contents rather than a path to the file. For reading, a StringIO will be used internally.

  • output_pdf (pathlib Path, string or None; default = None) – Path to the PDF file where the plot will be saved, if Path or string. (File is created.) If None, returns list of matplotlib figures.

  • which_beams (range or Sequence of int or strings) – Indices specifying which beams to plot. Order is the same as in the csv.

  • legends (None or list of strings, default = None) – Legends to be used for the various files.

Raises:

RuntimeError – If which_beams is not Sequence containing just integers or strings.

viperleed.calc.from_ase.rfactor_from_csv(beams_files, v0i, beams_file_is_content=(False, False), v0r_shift_range=(-3, 3), intpol_deg=5, intpol_step=0.5, return_beam_arrays=False)[source]

Compute the Pendry R-factor between two CSV files (in ViPErLEED format).

Read in two CSV files containing LEED-I(V) spectra and compute the mutual R-factor. The format of the files must be the same as accepted/generated by ViPErLEED.

Parameters:
  • beams_files (tuple of two string) – The two CSV files to be read in for R-factor calculation. These can contain either theoretical or experimental beams. If beams_file_is_content is set as True, the string is expected to contain the contents of the CSV file, otherwise it should be a path to the file.

  • v0i (float) – Imaginary part of the inner potential (in eV).

  • beams_file_is_content (tuple of two bool, default=(False, False)) – If either value is set to True, it is assumed that the corresponding element of beams_file contains a string with the CSV contents rather than a path to the file. For reading, a StringIO will be used internally.

  • v0r_shift_range (Sequence of two float, default=(-3, 3)) – During R-factor calculation, the two beam sets can be shifted in energy relative to one another. v0r_shift_range determines the minimum and maximum shifts (given in eV). The values should be integer multiples of intpol_step. If they are not, v0r_shift_range is expanded to the next integer multiples.

  • intpol_deg (int, default=5) – Either 3 or 5; degree of the interpolating spline polynomial

  • intpol_step (float, default=0.5) – Step size of the energy grid the beam data will be interpolated on (in eV).

  • return_beam_arrays (bool, default = False) – For debugging, you can optionally return the interpolated beams, y-functions, number of overlapping points & per-beam R-factors.

Returns:

  • best_R (float) – The R-factor between the two files. Returns the best R factor found after variation in v0r_shift_range.

  • best_v0r (float) – The v0r value for which the best R-factor was found.

viperleed.calc.from_ase.rot_mat_axis(axis, theta)[source]

Return a 3D rotation matrix by theta around axis.

Parameters:
  • axis (numpy.ndarray) – A vector parallel to the rotation axis. Shape (3,).

  • theta (float) – Rotation angle in degrees.

Returns:

rotation_matrix – Rotation matrix around axis by theta. The rotation is counter-clockwise when applied to column vectors from the left.

Return type:

numpy.ndarray

viperleed.calc.from_ase.rot_mat_x(theta)[source]

Return a rotation matrix around the x axis.

The rotation is positive, i.e., clockwise when looking along x, when applied from the left to column vectors. The same rotation for row vectors can be obtained by multiplying on the right with the transpose of the return value.

Parameters:

theta (float) – Angle of rotation in degrees.

Returns:

rot_mat – Shape (3, 3). Rotation matrix around the x axis.

Return type:

numpy.ndarray

viperleed.calc.from_ase.rot_mat_z(theta)[source]

Return a rotation matrix around the z axis.

The rotation is positive, i.e. clockwise when looking along z, when applied from the left to column vectors. The same rotation for row vectors can be obtained by multiplying on the right with the transpose of the return value.

Parameters:

theta (float) – Angle of rotation in degrees.

Returns:

rot_mat – Shape (3, 3). Rotation matrix around the z axis.

Return type:

numpy.ndarray

viperleed.calc.from_ase.run_from_ase(exec_path, ase_object, inputs_path=None, slab_transforms=(SlabTransform(orthogonal_matrix=None, uc_scaling=None, cut_cell_c_fraction=0.0),), cleanup_work=False)[source]

Perform a ViPErLEED calculation starting from an ase.Atoms object.

This function is a core part of the API for ViPErLEED. It allows to run a calculation from an ASE (Atomic Simulation Environment) Atoms object. The ase.Atoms object, which is expected to contain a simulation cell, will be cut at a specified height fraction, transformed into a ViPErLEED Slab and written to a POSCAR file. The ViPErLEED calculation is then executed on the structure using the input files found either in inputs_path or exec_path. The type of calculation performed depends on the value of the RUN parameter in the PARAMETERS file. If RUN = 1, this yields reference electron scattering amplitudes and intensities. These are collected in CSV files, which are returned as strings for further processing.

Parameters:
  • exec_path (string or Path) – Path where to execute the ViPErLEED calculation. The directory needs to exist. Input files will be copied from inputs_path, if provided. Otherwise it is assumed that all input files are already in exec_path. Temporary files will be stored in a subdirectory called work.

  • ase_object (ase.Atoms) – ASE-type object that contains the atoms for the slab to be used with ViPErLEED. This is expected to contain a simulation cell, which will be transformed using slab_transforms.

  • inputs_path (str or Path or None, optional) – Path to directory containing input files (PARAMETERS, VIBROCC, etc.) for ViPErLEED calculation. If provided, files will be copied from inputs_path to exec_path (silently overwriting existing files there!). If None or not given, no copying occurs. Default is None.

  • slab_transforms (SlabTransform or Sequence thereof, optional) – Sequence of transformation operations to be applied to the Slab constructed from ase_object. These include optional “change of basis” operations, optional uc_scaling for scaling the unit vector lengths, and a cut_cell_c_fraction for selecting only a portion of the ase_object. In addition, some special SlabTransform objects that swap unit cell axes can also be used (swap_a_b, swap_b_c, and swap_c_a). See help(SlabTransform) and help(swap_*_*) for further details. The transforms are applied in the order given. The cut_cell_c_fraction is taken only from the LAST SlabTransform given. Make sure that vectors a & b do not have any components in z after transformation as this is not allowed and will raise a ValueError. Default is a single SlabTransform(), corresponding to no transformation and no cut.

  • cleanup_work (bool, optional) – Whether the work directory created during execution of the calculation is to be removed at the end. Default is False.

Returns:

  • theobeams_file_str (str) – Contents of the file “THEOBEAMS.csv”, i.e. the scattering intensities as calculated in the reference calculation. An empty string is returned if no reference calculation was run.

  • amp_real_file_str (str) – Contents of the file “Complex_amplitudes_real.csv”, i.e., the real part of the scattering amplitudes as calculated in the reference calculation. An empty string is returned if no reference calculation was run or if the file was not found.

  • amp_imag_file_str (str) – Contents of the file “Complex_amplitudes_imag.csv”, i.e., the imaginary part of the scattering amplitudes as calculated in the reference calculation. An empty string is returned if no reference calculation was run or if the file was not found.

  • v0i (float) – Imaginary part of the inner potential as read from the PARAMETERS file. This is returned here since it is an input parameter for the R-factor calculation. The value is read from PARAMETERS and left unchanged - it is not a result of the calculation and is returned here for convenience only.

Raises:
  • FileNotFoundError – If exec_path is non-existent or not a directory.

  • FileNotFoundError – If no PARAMETERS file was found in exec_path, after potentially copying over from inputs_path

  • ValueError – If either of the first two unit cell vectors have a component perpendicular to the surface (i.e., along the z coordinate) after transformation of the unit cell

  • RuntimeError – If the ViPErLEED calculation fails.

Notes

The parameter slab_transform and its orthogonal_matrix attribute can be used to apply some useful transformations to the Slab. Examples include rotations around an axis and flipping the cell along c (in case the lower part has to be kept, rather than the upper one). A few special SlabTransform objects are also provided for swapping unit-cell vectors. Have a look at the matrices provided as part of this module (use as the orthogonal_matrix attribute of a SlabTransform): - rot_mat_x(theta) : Rotation matrix around x by theta (deg) - rot_mat_z(theta) : Rotation matrix around z by theta (deg) - rot_mat_axis(axis, theta): Rotation around axis by theta (deg) - flip_c_mat : Mirror matrix that flips the cell along c and the special SlabTransform objects: - swap_a_b : Swap vectors a & b (changes handedness) - swap_b_c : Swap vectors b & c (changes handedness) - swap_c_a : Swap vectors c & a (changes handedness) For applying multiple operations, the order does matter. You can: (i) Combine multiple orthogonal transformation matrices via matrix multiplication (@ symbol or np.dot) into a single SlabTransform object (into its orthogonal_matrix attribute): in this case the leftmost matrix is applied last. (ii) Provide a sequence of SlabTransform objects in slab_transforms: here transformations are applied in the order given (i.e., the leftmost is the applied first). (iii) A combination of the previous two.