Skip to content

setup #

Setup systems ready for calculations.

Functions:

  • setup_system

    Generate a set of molecule coordinate by using the PACKMOL package.

setup_system #

setup_system(
    components: list[tuple[str, int]],
    box_target_density: Quantity = 0.95 * _G_PER_ML,
    box_scale_factor: float = 1.1,
    box_padding: Quantity = 2.0 * angstrom,
    tolerance: Quantity = 2.0 * angstrom,
) -> tuple[Topology, Quantity]

Generate a set of molecule coordinate by using the PACKMOL package.

Parameters:

  • components (list[tuple[str, int]]) –

    A list of the form components[i] = (smiles_i, count_i) where smiles_i is the SMILES representation of component i and count_i is the number of corresponding instances of that component to create.

    If any SMILES is fully indexed an attempt will be made to create components with the same atom ordering as the indices.

  • box_target_density (Quantity, default: 0.95 * _G_PER_ML ) –

    Target mass density when approximating the box size for the final system with units compatible with g / mL.

  • box_scale_factor (float, default: 1.1 ) –

    The amount to scale the approximate box size by.

  • box_padding (Quantity, default: 2.0 * angstrom ) –

    The amount of extra padding to add to the box size to avoid PBC issues in units compatible with angstroms.

  • tolerance (Quantity, default: 2.0 * angstrom ) –

    The minimum spacing between molecules during packing in units compatible with angstroms.

Returns:

  • tuple[Topology, Quantity]

    A topology containing the molecules the coordinates were generated for and a unit [A] wrapped numpy array of coordinates with shape=(n_atoms, 3).

Source code in absolv/setup.py
def setup_system(
    components: list[tuple[str, int]],
    box_target_density: openmm.unit.Quantity = 0.95 * _G_PER_ML,
    box_scale_factor: float = 1.1,
    box_padding: openmm.unit.Quantity = 2.0 * openmm.unit.angstrom,
    tolerance: openmm.unit.Quantity = 2.0 * openmm.unit.angstrom,
) -> tuple[openff.toolkit.Topology, openmm.unit.Quantity]:
    """Generate a set of molecule coordinate by using the PACKMOL package.

    Args:
        components: A list of the form ``components[i] = (smiles_i, count_i)`` where
            ``smiles_i`` is the SMILES representation of component `i` and
            ``count_i`` is the number of corresponding instances of that component
            to create.

            If any SMILES is fully indexed an attempt will be made to create components
            with the same atom ordering as the indices.
        box_target_density: Target mass density when approximating the box size for the
            final system with units compatible with g / mL.
        box_scale_factor: The amount to scale the approximate box size by.
        box_padding: The amount of extra padding to add to the box size to avoid PBC
            issues in units compatible with angstroms.
        tolerance: The minimum spacing between molecules during packing in units
             compatible with angstroms.

    Returns:
        A topology containing the molecules the coordinates were generated for and
        a unit [A] wrapped numpy array of coordinates with shape=(n_atoms, 3).
    """

    packmol_path = shutil.which("packmol")

    if packmol_path is None:
        raise FileNotFoundError(errno.ENOENT, os.strerror(errno.ENOENT), "packmol")

    box_size = (
        _approximate_box_size_by_density(components, box_target_density)
        * box_scale_factor
    )

    molecules = {}

    for smiles, _ in components:
        if smiles in molecules:
            continue

        molecule = _molecule_from_smiles(smiles)
        molecule.generate_conformers(n_conformers=1)
        molecule.name = f"component-{len(molecules)}.xyz"
        molecules[smiles] = molecule

    with openff.utilities.temporary_cd():
        for molecule in molecules.values():
            molecule.to_file(molecule.name, "xyz")

        input_file_contents = _generate_input_file(
            [(molecules[smiles].name, count) for smiles, count in components],
            box_size,
            tolerance,
        )

        with open("input.txt", "w") as file:
            file.write(input_file_contents)

        with open("input.txt") as file:
            subprocess.run(packmol_path, stdin=file, check=True, capture_output=True)

        with open("output.xyz") as file:
            output_lines = file.read().splitlines(False)

    coordinates = (
        numpy.array(
            [
                [float(coordinate) for coordinate in coordinate_line.split()[1:]]
                for coordinate_line in output_lines[2:]
                if len(coordinate_line) > 0
            ]
        )
        * openmm.unit.angstrom
    )

    topology = openff.toolkit.Topology.from_molecules(
        [molecules[smiles] for smiles, count in components for _ in range(count)]
    )
    topology.box_vectors = numpy.eye(3) * (box_size + box_padding)

    return topology, coordinates