Skip to content

ase

Utility functions for running ASE calculators with ASE-based methods.

OptParams

Bases: TypedDict

Type hint for opt_params used throughout quacc.

fmax instance-attribute

fmax: float

fn_hook instance-attribute

fn_hook: Callable | None

max_steps instance-attribute

max_steps: int

optimizer class-attribute instance-attribute

optimizer: Optimizer = BFGS

optimizer_kwargs instance-attribute

optimizer_kwargs: OptimizerKwargs | None

run_kwargs instance-attribute

run_kwargs: dict[str, Any] | None

store_intermediate_results instance-attribute

store_intermediate_results: bool

OptimizerKwargs

Bases: TypedDict

Type hint for optimizer_kwargs in quacc.runners.ase.run_opt.

append_trajectory instance-attribute

append_trajectory: bool

restart instance-attribute

restart: Path | str | None

VibKwargs

Bases: TypedDict

Type hint for vib_kwargs in quacc.runners.ase.run_vib.

delta instance-attribute

delta: float

indices instance-attribute

indices: list[int] | None

nfree instance-attribute

nfree: int

run_calc

run_calc(
    atoms: Atoms,
    geom_file: str | None = None,
    copy_files: (
        SourceDirectory
        | dict[SourceDirectory, Filenames]
        | None
    ) = None,
    get_forces: bool = False,
) -> Atoms

Run a calculation in a scratch directory and copy the results back to the original directory. This can be useful if file I/O is slow in the working directory, so long as file transfer speeds are reasonable.

This is a wrapper around atoms.get_potential_energy(). Note: This function does not modify the atoms object in-place.

Parameters:

  • atoms (Atoms) –

    The Atoms object to run the calculation on.

  • geom_file (str | None, default: None ) –

    The filename of the log file that contains the output geometry, used to update the atoms object's positions and cell after a job. It is better to specify this rather than relying on ASE's atoms.get_potential_energy() function to update the positions, as this varies between codes.

  • copy_files (SourceDirectory | dict[SourceDirectory, Filenames] | None, default: None ) –

    Files to copy (and decompress) from source to the runtime directory.

  • get_forces (bool, default: False ) –

    Whether to use atoms.get_forces() instead of atoms.get_potential_energy().

Returns:

  • Atoms

    The updated Atoms object.

Source code in quacc/runners/ase.py
def run_calc(
    atoms: Atoms,
    geom_file: str | None = None,
    copy_files: SourceDirectory | dict[SourceDirectory, Filenames] | None = None,
    get_forces: bool = False,
) -> Atoms:
    """
    Run a calculation in a scratch directory and copy the results back to the original
    directory. This can be useful if file I/O is slow in the working directory, so long
    as file transfer speeds are reasonable.

    This is a wrapper around atoms.get_potential_energy(). Note: This function
    does not modify the atoms object in-place.

    Parameters
    ----------
    atoms
        The Atoms object to run the calculation on.
    geom_file
        The filename of the log file that contains the output geometry, used to
        update the atoms object's positions and cell after a job. It is better
        to specify this rather than relying on ASE's
        atoms.get_potential_energy() function to update the positions, as this
        varies between codes.
    copy_files
        Files to copy (and decompress) from source to the runtime directory.
    get_forces
        Whether to use `atoms.get_forces()` instead of `atoms.get_potential_energy()`.

    Returns
    -------
    Atoms
        The updated Atoms object.
    """
    # Copy atoms so we don't modify it in-place
    atoms = copy_atoms(atoms)

    # Perform staging operations
    tmpdir, job_results_dir = calc_setup(atoms, copy_files=copy_files)

    # Run calculation
    try:
        if get_forces:
            atoms.get_forces()
        else:
            atoms.get_potential_energy()
    except Exception as exception:
        terminate(tmpdir, exception)

    # Most ASE calculators do not update the atoms object in-place with a call
    # to .get_potential_energy(), which is important if an internal optimizer is
    # used. This section is done to ensure that the atoms object is updated to
    # the final geometry if `geom_file` is provided.
    # Note: We have to be careful to make sure we don't lose the calculator
    # object, as this contains important information such as the parameters
    # and output properties (e.g. final magnetic moments).
    if geom_file:
        atoms_new = read(zpath(tmpdir / geom_file))
        if isinstance(atoms_new, list):
            atoms_new = atoms_new[-1]

        # Make sure the atom indices didn't get updated somehow (sanity check).
        # If this happens, there is a serious problem.
        if (
            np.array_equal(atoms_new.get_atomic_numbers(), atoms.get_atomic_numbers())
            is False
        ):
            raise ValueError("Atomic numbers do not match between atoms and geom_file.")

        atoms.positions = atoms_new.positions
        atoms.cell = atoms_new.cell

    # Perform cleanup operations
    calc_cleanup(atoms, tmpdir, job_results_dir)

    return atoms

run_opt

run_opt(
    atoms: Atoms,
    relax_cell: bool = False,
    fmax: float = 0.01,
    max_steps: int = 1000,
    optimizer: Optimizer = BFGS,
    optimizer_kwargs: OptimizerKwargs | None = None,
    store_intermediate_results: bool = False,
    fn_hook: Callable | None = None,
    run_kwargs: dict[str, Any] | None = None,
    copy_files: (
        SourceDirectory
        | dict[SourceDirectory, Filenames]
        | None
    ) = None,
) -> Optimizer

Run an ASE-based optimization in a scratch directory and copy the results back to the original directory. This can be useful if file I/O is slow in the working directory, so long as file transfer speeds are reasonable.

This is a wrapper around the optimizers in ASE. Note: This function does not modify the atoms object in-place.

Parameters:

  • atoms (Atoms) –

    The Atoms object to run the calculation on.

  • relax_cell (bool, default: False ) –

    Whether to relax the unit cell shape and volume.

  • fmax (float, default: 0.01 ) –

    Tolerance for the force convergence (in eV/A).

  • max_steps (int, default: 1000 ) –

    Maximum number of steps to take.

  • optimizer (Optimizer, default: BFGS ) –

    Optimizer class to use.

  • optimizer_kwargs (OptimizerKwargs | None, default: None ) –

    Dictionary of kwargs for the optimizer. Takes all valid kwargs for ASE Optimizer classes. Refer to _set_sella_kwargs for Sella-related kwargs and how they are set.

  • store_intermediate_results (bool, default: False ) –

    Whether to store the files generated at each intermediate step in the optimization. If enabled, they will be stored in a directory named stepN where N is the step number, starting at 0.

  • fn_hook (Callable | None, default: None ) –

    A custom function to call after each step of the optimization. The function must take the instantiated dynamics class as its only argument.

  • run_kwargs (dict[str, Any] | None, default: None ) –

    Dictionary of kwargs for the run() method of the optimizer.

  • copy_files (SourceDirectory | dict[SourceDirectory, Filenames] | None, default: None ) –

    Files to copy (and decompress) from source to the runtime directory.

Returns:

  • Optimizer

    The ASE Optimizer object.

Source code in quacc/runners/ase.py
def run_opt(
    atoms: Atoms,
    relax_cell: bool = False,
    fmax: float = 0.01,
    max_steps: int = 1000,
    optimizer: Optimizer = BFGS,
    optimizer_kwargs: OptimizerKwargs | None = None,
    store_intermediate_results: bool = False,
    fn_hook: Callable | None = None,
    run_kwargs: dict[str, Any] | None = None,
    copy_files: SourceDirectory | dict[SourceDirectory, Filenames] | None = None,
) -> Optimizer:
    """
    Run an ASE-based optimization in a scratch directory and copy the results back to
    the original directory. This can be useful if file I/O is slow in the working
    directory, so long as file transfer speeds are reasonable.

    This is a wrapper around the optimizers in ASE. Note: This function does not
    modify the atoms object in-place.

    Parameters
    ----------
    atoms
        The Atoms object to run the calculation on.
    relax_cell
        Whether to relax the unit cell shape and volume.
    fmax
        Tolerance for the force convergence (in eV/A).
    max_steps
        Maximum number of steps to take.
    optimizer
        Optimizer class to use.
    optimizer_kwargs
        Dictionary of kwargs for the optimizer. Takes all valid kwargs for ASE
        Optimizer classes. Refer to `_set_sella_kwargs` for Sella-related
        kwargs and how they are set.
    store_intermediate_results
        Whether to store the files generated at each intermediate step in the
        optimization. If enabled, they will be stored in a directory named
        `stepN` where `N` is the step number, starting at 0.
    fn_hook
        A custom function to call after each step of the optimization.
        The function must take the instantiated dynamics class as
        its only argument.
    run_kwargs
        Dictionary of kwargs for the run() method of the optimizer.
    copy_files
        Files to copy (and decompress) from source to the runtime directory.

    Returns
    -------
    Optimizer
        The ASE Optimizer object.
    """
    # Copy atoms so we don't modify it in-place
    atoms = copy_atoms(atoms)

    # Perform staging operations
    tmpdir, job_results_dir = calc_setup(atoms, copy_files=copy_files)

    # Set defaults
    optimizer_kwargs = recursive_dict_merge(
        {
            "logfile": "-" if SETTINGS.DEBUG else tmpdir / "opt.log",
            "restart": tmpdir / "opt.json",
        },
        optimizer_kwargs,
    )
    run_kwargs = run_kwargs or {}

    # Check if trajectory kwarg is specified
    if "trajectory" in optimizer_kwargs:
        msg = "Quacc does not support setting the `trajectory` kwarg."
        raise ValueError(msg)

    # Handle optimizer kwargs
    if optimizer.__name__.startswith("SciPy"):
        optimizer_kwargs.pop("restart", None)
    elif optimizer.__name__ == "Sella":
        _set_sella_kwargs(atoms, optimizer_kwargs)
    elif optimizer.__name__ == "IRC":
        optimizer_kwargs.pop("restart", None)

    # Define the Trajectory object
    traj_file = tmpdir / "opt.traj"
    traj = Trajectory(traj_file, "w", atoms=atoms)
    optimizer_kwargs["trajectory"] = traj

    # Set volume relaxation constraints, if relevant
    if relax_cell and atoms.pbc.any():
        atoms = FrechetCellFilter(atoms)

    # Run optimization
    try:
        with traj, optimizer(atoms, **optimizer_kwargs) as dyn:
            if optimizer.__name__.startswith("SciPy"):
                # https://gitlab.com/ase/ase/-/issues/1475
                dyn.run(fmax=fmax, steps=max_steps, **run_kwargs)
            else:
                for i, _ in enumerate(
                    dyn.irun(fmax=fmax, steps=max_steps, **run_kwargs)
                ):
                    if store_intermediate_results:
                        _copy_intermediate_files(
                            tmpdir,
                            i,
                            files_to_ignore=[
                                traj_file,
                                optimizer_kwargs["restart"],
                                optimizer_kwargs["logfile"],
                            ],
                        )
                    if fn_hook:
                        fn_hook(dyn)
    except Exception as exception:
        terminate(tmpdir, exception)

    # Store the trajectory atoms
    dyn.traj_atoms = read(traj_file, index=":")

    # Perform cleanup operations
    calc_cleanup(get_final_atoms_from_dynamics(dyn), tmpdir, job_results_dir)

    return dyn

run_vib

run_vib(
    atoms: Atoms,
    vib_kwargs: VibKwargs | None = None,
    copy_files: (
        SourceDirectory
        | dict[SourceDirectory, Filenames]
        | None
    ) = None,
) -> Vibrations

Run an ASE-based vibration analysis in a scratch directory and copy the results back to the original directory. This can be useful if file I/O is slow in the working directory, so long as file transfer speeds are reasonable.

This is a wrapper around the vibrations module in ASE. Note: This function does not modify the atoms object in-place.

Parameters:

Returns:

Source code in quacc/runners/ase.py
def run_vib(
    atoms: Atoms,
    vib_kwargs: VibKwargs | None = None,
    copy_files: SourceDirectory | dict[SourceDirectory, Filenames] | None = None,
) -> Vibrations:
    """
    Run an ASE-based vibration analysis in a scratch directory and copy the results back
    to the original directory. This can be useful if file I/O is slow in the working
    directory, so long as file transfer speeds are reasonable.

    This is a wrapper around the vibrations module in ASE. Note: This function
    does not modify the atoms object in-place.

    Parameters
    ----------
    atoms
        The Atoms object to run the calculation on.
    vib_kwargs
        Dictionary of kwargs for the [ase.vibrations.Vibrations][] class.
    copy_files
        Files to copy (and decompress) from source to the runtime directory.

    Returns
    -------
    Vibrations
        The updated Vibrations module
    """
    # Copy atoms so we don't modify it in-place
    atoms = copy_atoms(atoms)

    # Set defaults
    vib_kwargs = vib_kwargs or {}

    # Perform staging operations
    tmpdir, job_results_dir = calc_setup(atoms, copy_files=copy_files)

    # Run calculation
    vib = Vibrations(atoms, name=str(tmpdir / "vib"), **vib_kwargs)
    try:
        vib.run()
    except Exception as exception:
        terminate(tmpdir, exception)

    # Summarize run
    vib.summary(log=sys.stdout if SETTINGS.DEBUG else str(tmpdir / "vib_summary.log"))

    # Perform cleanup operations
    calc_cleanup(vib.atoms, tmpdir, job_results_dir)

    return vib