wfl.generate package#

Submodules#

wfl.generate.atoms_and_dimers module#

wfl.generate.atoms_and_dimers.isolated_atom_from_e0(outputs, e0_dict, cell_size, energy_key='energy', extra_info=None)#

Write single atoms with energy from e0 dictionary

Parameters
  • outputs (Configset_out) –

  • e0_dict (dict) –

  • cell_size (float) – unit cell size, should be greater than cutoff

  • energy_key (dict_key, default 'energy') –

  • extra_info (dict) – extra dict to fill into atoms.info

wfl.generate.atoms_and_dimers.prepare(outputs, atomic_numbers, bond_lengths=None, dimer_n_steps=41, dimer_factor_range=(0.65, 2.5), dimer_box_scale=6.0, extra_info=None, do_isolated_atoms=True, max_cutoff=None, fixed_cell=False)#

Prepare dimer and single atom configurations

Parameters
  • outputs (OutputSpec) – target for atom and dimer configurations

  • atomic_numbers (list(int)) – list of atomic numbers

  • bond_lengths (dict, None, default={Any -> 1.0}) – bond lengths by species, default is to assign 1.0 to all elements

  • dimer_n_steps (int) – number of dimer steps to take

  • dimer_factor_range (tuple, default=(0.65, 2.5)) – Range of factors of bond length to use for dimer separations. If bond_length is not given, this becomes the exact range of distance for all pairs.

  • dimer_box_scale (float, default=6.0) – dimer box size, as scale of bond length

  • extra_info (dict) – extra dict to fill into atoms.info

  • do_isolated_atoms (bool, default=True) – create isolated atoms in boxes as well, cell is bond_len * dimer_box_scale in all three directions or fixed cell requires max_cutoff

  • max_cutoff (float, default None) – max cutoff, needed to ensure that isolate_atom boxes are large enough for gap_fit to accept

  • fixed_cell (ase.cell.Cell, list, bool, default=False) – switch to use a fixed cell, use the cell given. Use anything that ase.cell.Cell understands

Return type

ConfigSet in with configs

wfl.generate.buildcell module#

wfl.generate.buildcell.buildcell(*args, **kwargs)#

Creates atomic configurations by repeatedly running buildcell

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • buildcell_cmdbuildcell executable, including path if needed

  • buildcell_input (str) – input file contents to buildcell, in the form of a single string

  • extra_info (dict, default {}) – extra fields to place into created atoms info dict

  • perturbation (float, default 0.0) – magnitude of perturbation to atomic positions with Atoms.rattle()

  • skip_failures (bool, default True) – skip failed buildcell calls

  • symprec (float, default 0.01) – precision for symmetry check

  • verbose (bool, default False) – print some running info

Return type

list of Atoms created by buildcell

autopara_infoAutoParaInfo / dict, optional

information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.buildcell.conv_buildcell_out(buildcell_output)#

Convert from buildcell output to atoms object(s)

Parameters

buildcell_output (str) – stdout of buildcell run

Returns

ats – list of created Atoms objects

Return type

list(Atoms)

wfl.generate.buildcell.create_input(z, vol_per_atom, bond_lengths, composition=None, RSS_min_vol_factor=0.5, vol_range=(0.95, 1.05), min_sep_factor=0.9, symmops='1-8', natom=(6, 24), odd=None, verbose=False)#

Create an input string suitable for a random structure search (RSS) via AIRSS’ buildcell.

Parameters
  • z (int / list(int)) – Atomic number(s) to be present in the system.

  • vol_per_atom (float / list(float)) – Volume to be reserved per atom in the system (thus, affects overall unit cell volume/system density).

  • bond_lengths (float / list(float)) – Defines the separation between atoms.

  • composition (list(int), default None) – Defines the stoichiometric ratio between species, in case multiple species are present (i.e. z is a list)

  • RSS_min_vol_factor (float, default 0.5) – Together with vol_per_atom defines the minimum volume to be reserved per atom.

  • vol_range (tuple(2), default (0.95, 1.05)) – Defines the variations allowed in the overall unit cell volume.

  • min_sep_factor (float, default 0.9) – Together with bond_lengths defines the minimum distance between atoms.

  • symmops (str, default '1-8') – Build structures have a specified number of symmetry operations. For crystals, the allowed values are (1,2,3,4,6,8,12,16,24,48). For clusters (indicated with #CLUSTER), the allowed values are (1,2,3,5,4,6,7,8,9,10,11,12,24). Ranges are allowed (see default).

  • natom (tuple(2), default (6, 24)) – Defines the range of total numbers of atoms to be present in build structures.

  • odd (str, default None) – Control odd/even numbers of atoms in cell. Use "only" to only produce cells with odd numbers of atoms, "also" to produce both odd and even, and None for even only.

  • verbose (bool, default False) – verbose output

Returns

output

Return type

str suitable as input file contents to buildcell

wfl.generate.minimahopping module#

wfl.generate.minimahopping.minimahopping(*args, **kwargs)#

runs a structure optimization

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • calculator (Calculator / (initializer, args, kwargs)) – ASE calculator or routine to call to create calculator

  • Ediff0 (float, default 1 (eV)) – initial energy acceptance threshold

  • T0 (float, default 1000 (K)) – initial MD temperature

  • minima_threshold (float, default 0.5 (A)) – threshold for identical configs

  • mdmin (int, default 2) – criteria to stop MD simulation (number of minima)

  • fmax (float, default 1 (eV/A)) – max force for optimizations

  • timestep (float, default 1 (fs)) – timestep for MD simulations

  • totalsteps (int, default 10) – number of steps

  • skip_failures (bool, default True) – just skip optimizations that raise an exception

  • workdir (str/Path default None) – workdir for saving files

  • opt_kwargs – keyword arguments for MinimaHopping

  • rng (numpy.random.Generator, default None) – random number generator to use (needed for pressure sampling, initial temperature, or Langevin dynamics)

  • _autopara_per_item_info (dict) – INTERNALLY used by autoparallelization framework to make runs reproducible (see wfl.autoparallelize.autoparallelize() docs)

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.neb module#

wfl.generate.neb.NEB(*args, **kwargs)#

runs a structure optimization

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • calculator (Calculator / (initializer, args, kwargs)) – ASE calculator or initializer and arguments to call to create calculator

  • fmax (float, default 5e-2) – force convergence tolerance

  • steps (int, default 1000) – max number of steps

  • traj_step_interval (int, default 1) – if present, interval between trajectory snapshots

  • traj_subselect ("last_converged", default None) – rule for sub-selecting configs from the full trajectory. Currently implemented: “last_converged”, which takes the last config, if converged.

  • skip_failures (bool, default True) – just skip optimizations that raise an exception

  • verbose (bool, default False) – verbose output optimisation logs are not printed unless this is True

  • update_config_type (bool, default True) – append at.info[‘optimize_config_type’] at.info[‘config_type’]

  • neb_kwargs – keyword arguments for DyNEB and FIRE

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.neb.subselect_from_traj(traj, subselect=None)#

Sub-selects configurations from trajectory.

Parameters

subselect (int or string, default None) –

  • None: full trajectory is returned

  • int: (not implemented) how many samples to take from the trajectory.

  • str: specific method

    • ”last_converged”: returns [last_config] if converged, or None if not.

wfl.generate.normal_modes module#

Generate configs with normal mode sampling, generally useful for small molecules to learn around their relaxed state.

class wfl.generate.normal_modes.NormalModes(nm_atoms, prop_prefix)#

Bases: object

derive_normal_mode_info(calculator, parallel_hessian=True)#

Get normal mode information using numerical hessian

Parameters
  • calculator (Calculator / (initializer, args, kwargs)) – ASE calculator or routine to call to create calculator

  • parallel_hessian (bool, default=True) – whether to parallelize 6N calculations needed for approximating the Hessian.

static evals_to_freqs(evals)#

Converts from eigenvalues in eV/Å^2/amu to frequencies (square root of eigenvalue) in eV. Negative frequencies are shorthand for imaginary frequencies.

static evecs_to_modes(evecs, masses=None, inverse_m=None)#

converts from mass-weighted 3N-long eigenvector to 3N displacements in Cartesian coordinates

static freqs_to_evals(freqs)#

Converts from frequencies (sqrt of eigenvalue) in eV to eigenvalues in eV/Å^2/amu. Negative frequencies are shorthand for imaginary frequencies.

static modes_to_evecs(modes, masses=None, inverse_m=None)#

converts 3xN cartesian displacements to 1x3N mass-weighted eigenvectors

num_hess_delta = 0.01#
sample_normal_modes(sample_size, temp=None, energies_for_modes=None, normal_mode_numbers='all', info_to_keep='default', arrays_to_keep=None)#

Randomly displace equilibrium structure’s atoms along given normal modes so that normal mode energies follow Boltzmann distribution at a given temperature.

Notes

Coefficients for scaling the individual mass-weighted normal coordinates are sampled from a normal distribution with standard deviation of eigenvalue / (kB * temp). The harmonic energy distribution then follows a generalised gamma distribution with degrees of freedom = len(normal_mode_numbers) or shape = len(normal_mode_numbers) / 2.

Returned atoms have an atoms.info[‘nm_energy’] entry with the corresponding harmonic energy

Parameters
  • sample_size (int) – How many randomly perturbed structures to return

  • temp (float, default None) – Temperature for the displacement (putting kT energy into each mode) alternative to energies_for_modes

  • energies_for_modes (list(float), default: None) – list of energies (e.g. kT) for each of the normal modes for generating displacement magnitudes. Alternative to temp. Length either 3N - 6 if normal_mode_numbers == ‘all’ or matches len(normal_mode_numbers).

  • normal_mode_numbers (int / list(int) / 'all', default "all") – List of normal mode numbers to displace along. Alternatively if “all” is selected, all but first six (rotations and translations) normal modes are used.

  • info_to_keep (list(str) / None / 'default', default 'default') – list of Atoms.info.keys() to keep, defaults to ‘config_type’ only

  • arrays_to_keep (list(str), default None) – list of Atoms.arrays.keys() entries to keep

Return type

list(Atoms) of randomly perturbed structures

summary()#

Prints all vibrational frequencies.

view(prefix='nm', output_dir='normal_modes', normal_mode_numbers='all', temp=300, nimages=16)#

writes out xyz files with oscillations along each of the normal modes

Parameters
  • prefix (str, default "nm") – Prefix normal mode files

  • output_dir (str, default "normal_modes") – Directory for outputs

  • normal_mode_numbers (str / list(int) / np.array(int), default "all") – List of normal mode numbers to write out or “all” to save all

  • temp (float, default 300) – Temperature for the oscillation

  • nimages (int, default 32) – Number of structures per oscillation/output file

wfl.generate.normal_modes.generate_normal_modes_parallel_atoms(*args, **kwargs)#

Get normal mode information for all atoms in the input

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • calculator (Calculator / (initializer, args, kwargs)) – ASE calculator or routine to call to create calculator

  • prop_prefix (str / None) – prefix for normal_mode_frequencies and normal_mode_displacements stored in atoms.info/arrays

  • parallel_hessian (bool, default=True) – whether to parallelize 6N calculations needed for Hessian approximation.

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.normal_modes.generate_normal_modes_parallel_hessian(inputs, outputs, calculator, prop_prefix)#
wfl.generate.normal_modes.sample_normal_modes(inputs, outputs, temp, sample_size, prop_prefix, info_to_keep=None, arrays_to_keep=None)#

Multiple times displace along normal modes for all atoms in input

Parameters
  • inputs (Atoms / list(Atoms) / ConfigSet) – Structures with normal mode information (eigenvalues & eigenvectors)

  • outputs (OutputSpec) –

  • temp (float) – Temperature for normal mode displacements

  • sample_size (int) – Now many perturbed structures per input structure to return

  • prop_prefix (str / None) – prefix for normal_mode_frequencies and normal_mode_displacements stored in atoms.info/arrays

  • info_to_keep (str, default "config_type") – string of Atoms.info.keys() to keep

  • arrays_to_keep (str, default None) – string of Atoms.arrays.keys() entries to keep

wfl.generate.optimize module#

wfl.generate.optimize.optimize(*args, **kwargs)#

runs a structure optimization

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • calculator (Calculator / (initializer, args, kwargs)) – ASE calculator or routine to call to create calculator

  • fmax (float, default 1e-3) – force convergence tolerance

  • smax (float, default None) – stress convergence tolerance, default from fmax

  • steps (int, default 1000) – max number of steps

  • pressure (None / float / tuple) – applied pressure distribution (GPa), as parsed by wfl.utils.pressure.sample_pressure()

  • stress_mask (None / list(bool)) – mask for stress components to pass to variable-cell filter

  • keep_symmetry (bool, default True) – constrain symmetry to maintain initial

  • traj_step_interval (int, default 1) – if present, interval between trajectory snapshots

  • traj_subselect ("last_converged", default None) – rule for sub-selecting configs from the full trajectory. Currently implemented: “last_converged”, which takes the last config, if converged.

  • skip_failures (bool, default True) – just skip optimizations that raise an exception

  • verbose (bool, default False) – verbose output optimisation logs are not printed unless this is True

  • update_config_type (bool, default True) – append at.info[‘optimize_config_type’] at.info[‘config_type’]

  • opt_kwargs – keyword arguments for PreconLBFGS

  • rng (numpy.random.Generator, default None) – random number generator to use (needed for pressure sampling, initial temperature, or Langevin dynamics)

  • _autopara_per_item_info (dict) – INTERNALLY used by autoparallelization framework to make runs reproducible (see wfl.autoparallelize.autoparallelize() docs)

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.optimize.subselect_from_traj(traj, subselect=None)#

Sub-selects configurations from trajectory.

Parameters

subselect (int or string, default None) –

  • None: full trajectory is returned

  • int: (not implemented) how many samples to take from the trajectory.

  • str: specific method

    • ”last”: returns [last_config]

    • ”last_converged”: returns [last_config] if converged, or None if not.

wfl.generate.phonopy module#

wfl.generate.phonopy.phonopy(*args, **kwargs)#

create displaced configs with phonopy or phono3py for each structure in inputs

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • displacements (list(float)) – list of displacement magnitudes to use

  • strain_displs (list(float)) – list of deformation gradient magnitudes to use

  • ph2_supercell (int or int 3-vector or 3x3 matrix, default None) – supercell to use for harmonic force constant displacements

  • ph2_supercell – if not None, supercell to use for cubic force constant displacements

  • pair_cutoff (float, default None) – if not None, max distance between pairs of atoms displaced for cubic force constants

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.smiles module#

wfl.generate.smiles.smi_to_atoms(smi, useBasicKnowledge=True, useExpTorsionAnglePrefs=True, randomSeed=-1)#

Converts smiles to 3D Atoms object

wfl.generate.smiles.smiles(*args, **kwargs)#

Creates atomic configurations by repeatedly running smi_to_xyz, I/O with OutputSpec.

Parameters
  • inputs (iterable(SMILES string)) – input quantities of type SMILES string

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • useBasicKnowledge (bool, default True) – impose basic knowledge such as flat aromatic rings

  • useExpTorsionAnglePrefs (bool, default True) – impose experimental torsion angle preferences

  • extra_info (dict, default {}) – extra fields to place into created atoms info dict

  • randomSeed (int, default -1) – RDKit EmbedMolecule random seed for reproducibility

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.supercells module#

wfl.generate.supercells.antisite(*args, **kwargs)#

make antisites in largest bulk-like supercells

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • max_n_atoms (int) – maximum number of atoms allowed in a supercell

  • pert (float, default 0.0) – magnitude of random displacement

  • primitive (bool, default True) – reduce input Atoms to primitive cell using spglib

  • symprec (float, default 1.0e-3) – symprec for primitive lattice check

  • n_antisite (int, default 1) – number of antisites to create in each supercell

  • cluster_r (float, default 0.0) – if > 0.0, multiply by 1st nn distance of initial antisite, and make antisites in cluster of atoms within this distance of initial antisite.

  • Zs (iterable(int), default None) – list of atomic numbers to be used to create antisites. If none or length 0, use set of species already present in each config for itself.

  • rng (numpy.random.Generator) – random number generator

  • _autopara_per_item_info (list) – internal use

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.supercells.interstitial(*args, **kwargs)#

make interstitials in largest bulk-like supercells

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • max_n_atoms (int) – maximum number of atoms allowed in a supercell

  • pert (float, default 0.0) – magnitude of random displacement

  • interstitial_probability_radius_exponent (float, default 3.0) – probability of selecting an interstitial void is proportional to its radius to this power

  • primitive (bool, default True) – reduce input Atoms to primitive cell using spglib

  • symprec (float, default 1.0e-3) – symprec for primitive lattice check

  • rng (numpy.random.Generator) – random number generator

  • _autopara_per_item_info (list) – internal use

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.supercells.largest_bulk(*args, **kwargs)#

make largest bulk-like supercells

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • max_n_atoms (int) – maximum number of atoms allowed in a supercell

  • pert (float, default 0.0) – magnitude of random displacement

  • primitive (bool, default True) – reduce input Atoms to primitive cell using spglib

  • symprec (float, default 1.0e-3) – symprec for primitive lattice check

  • ase_optimal (bool, default False) – use ase.build.supercells.find_optimal_cell_shape

  • rng (numpy.random.Generator) – random number generator

  • _autopara_per_item_info (list) – internal use

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.supercells.surface(*args, **kwargs)#

make surface supercells

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • max_n_atoms (int) – maximum number of atoms allowed in a supercell

  • min_thickness (float) – duplicate slab in normal direction to exceed this minimum thickness

  • vacuum (float) – thickness of vacuum layer to add

  • simple_cut (bool, default False) – if true, surface is formed by two primitive lattice vectors, otherwise use two random vectors with indices up to +/- max_surface_cell_indices

  • max_surface_cell_indices (int, default 1) – maximum indices used to make surfaces when simple_cut=False

  • duplicate_in_plane (bool, default True) – duplicate surface cell in-plane if allowed by max_n_atoms

  • pert (float, default 0.0) – magnitude of random displacement

  • primitive (bool, default True) – reduce input Atoms to primitive cell using spglib

  • symprec (float, default 1.0e-3) – symprec for primitive lattice check

  • rng (numpy.random.Generator) – random number generator

  • _autopara_per_item_info (list) – internal use

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.supercells.vacancy(*args, **kwargs)#

make vacancies in largest bulk-like supercells

Parameters
  • inputs (iterable(Atoms)) – input quantities of type Atoms

  • outputs (OutputSpec or None) – where to write output atomic configs, or None for no output (i.e. only side-effects)

  • max_n_atoms (int) – maximum number of atoms allowed in a supercell

  • pert (float, default 0.0) – magnitude of random displacement

  • primitive (bool, default True) – reduce input Atoms to primitive cell using spglib

  • symprec (float, default 1.0e-3) – symprec for primitive lattice check

  • n_vac (int, default 1) – number of vacancies to create in each supercell

  • cluster_r (float, default 0.0) – if > 0.0, multiply by 1st nn distance of initial vacancy, and make vacancies in cluster of atoms within this distance of initial vacancy.

  • rng (numpy.random.Generator) – random number generator

  • _autopara_per_item_info (list) – internal use

  • autopara_info (AutoParaInfo / dict, optional) – information for automatic parallelization

Returns

co – output configs

Return type

ConfigSet

wfl.generate.utils module#

wfl.generate.utils.config_type_append(at, config_type)#

Module contents#