API reference

The AtomPacker is a Python package for packing nanoclusters into supramolecular cages. It provides a set of tools for detecting cavities in molecular structures, packing nanoclusters into these cavities, and exporting the resulting structures.

class AtomPacker.Cage[source]

A class representing a supramolecular cage loaded from a structure file.

The AtomPacker.Cage class is used to store the atoms of a macromolecular structure.

property atomic: ndarray | None

Get the coordinates of the atoms in the molecular structure.

Returns:

An array containing the atomic information. The array has the following columns: resnum, chain, resname, element, x, y, z, radius. If the cage is not loaded, returns None.

Return type:

numpy.ndarray | None

property centroid: ndarray | None

Get the centroid of the cage structure.

Returns:

An array containing the xyz coordinates of the centroid of the cage structure. If the cage is not loaded, returns None.

Return type:

numpy.ndarray | None

Raises:

ValueError – If the cage is not loaded.

property coordinates: ndarray | None

Get the coordinates of the atoms in the molecular structure.

Returns:

An array containing the atomic coordinates. The array has the following columns: x, y, z. If the cage is not loaded, returns None.

Return type:

numpy.ndarray | None

detect_cavity(step: float = 0.6, probe_in: float = 1.4, probe_out: float = 10.0, removal_distance: float = 2.4, volume_cutoff: float = 100.0, surface: str = 'SES', nthreads: int | None = None, verbose: bool = False, **kwargs: dict[str, object]) None[source]

Detect the cavity in the macromolecular structure.

Parameters:
  • step (float, optional) – The grid spacing for the cavity detection algorithm.

  • probe_in (float, optional) – The radius of the inner probe.

  • probe_out (float, optional) – The radius of the outer probe.

  • removal_distance (float, optional) – The distance to remove atoms from the cavity.

  • volume_cutoff (float, optional) – The volume cutoff for the cavity.

  • surface (str, optional) – Surface representation. Keywords options are SES (Solvent Excluded Surface) or SAS (Solvent Accessible Surface), by default SES.

  • nthreads (int, optional) – Number of threads, by default None. If None, the number of threads is os.cpu_count() - 1. If -1, uses all available threads.

  • verbose (bool, optional) – Print extra information to standard output, by default False.

  • kwargs (dict[str, object], optional) – Additional keyword arguments to pass to the cavity detection algorithm (pyKVFinder.detect) of pyKVFinder package.

Returns:

The cavity object. The cavity object contains the grid points of the cavity, the step size, the Probe In radius, the Probe Out radius, the removal distance, the volume cutoff, and the vertices (origin, X-axis, Y-axis, Z-axis) of the grid. Cavity points in the 3D grid (grid[nx][ny][nz]). Cavities array has integer labels in each position, that are:

  • -1: solvent points;

  • 0: cage points;

  • 1: empty space points;

  • >=2: cavity points.

The empty space points are regions that do not meet the chosen volume cutoff to be considered a cavity.

Return type:

Cavity

Raises:

ValueError – If the cage is not loaded.

load(filename: str, vdw: dict[str, float] | None = None) None[source]

Load a supramolecular cage structure file into the MDAnalysis .Univese object.

Parameters:
  • filename (str) – The filename of the structure file. The file format is determined by the suffix. Supported formats are: .cif, .mol2, .pdb, .xyz.

  • vdw (dict[str, float] | None, optional) – A dictionary containing the van der Waals radii for each atom type, by default None. If None, the van der Waals radii are looked up from the pyKVFinder package.

Returns:

universe – The MDAnalysis universe object containing the structure.

Return type:

MDAnalysis.Univese

Raises:

ValueError – If the file format is not supported.

pack(atom_type: str, lattice_type: str, a: float | None = None, b: float | None = None, c: float | None = None, clashing_tolerance: float = 0.0, angles: ndarray | list[float] = [0.0], translations: ndarray | list[float] = [0.0], optsave: bool = False, optdir: str | None = None, nthreads: int | None = None, verbose: bool = False) None[source]

Pack a cluster of atoms into the cage structure.

Parameters:
  • atom_type (str) – Atom type for the cluster.

  • lattice_type (str) – Lattice type for the cluster. Supported types: ‘bcc’, ‘fcc’, ‘hcp’, ‘sc’.

  • a (float | None, optional) – Lattice constant ‘a’. If None, uses values from AtomPacker or ASE.

  • b (float | None, optional) – Lattice constant ‘b’. If None, uses values from AtomPacker or ASE.

  • c (float | None, optional) – Lattice constant ‘c’. If None, uses values from AtomPacker or ASE.

  • clashing_tolerance (float, optional) – Minimum allowed distance (Å) between cluster and cage atoms. Default is 0.0.

  • angles (numpy.ndarray | list[float], optional) – Rotation angles for cluster optimization. If not specified, no optimization is performed. Default is [0.0]. If specified, angles should be a list or numpy array of angles in degrees. Example: [-75, -50, -25, 0, 25, 50, 75].

  • translations (numpy.ndarray | list[float], optional) – Translation values for cluster optimization. If not specified, no optimization is performed. Default is [0.0]. If specified, translations should be a list or numpy array of translation values in Angstroms. Example: [-0.2, 0.0, 0.2].

  • optsave (bool, optional) – If True, saves each optimization step as a PDB file. Default is False.

  • optdir (str | None, optional) – Directory to save files. If None, uses current working directory.

  • nthreads (int | None, optional) – Number of threads to use for parallel processing. If None, uses os.cpu_count() - 1. Default is None. If -1, uses all available threads.

  • verbose (bool, optional) – If True, prints detailed information during processing (default is False).

Raises:

ValueError – If the cage is not loaded, cavity is not detected, or clashing_tolerance < 0.

preview(show_cavity: bool = False, show_cluster: bool = False, show_openings: bool = False, renderer: str = 'browser', **kwargs: dict[str, object]) None[source]

Preview the cage system (cage, cavity, and cluster) in a 3D viewer.

Parameters:
  • show_cavity (bool, optional) – Show the cavity in the 3D viewer, by default False.

  • show_cluster (bool, optional) – Show the cluster in the 3D viewer, by default False.

  • show_openings (bool, optional) – Show the openings in the cavity in the 3D viewer, by default False.

  • renderer (str, optional) – The renderer to use for the 3D viewer. Supported renderers are ‘browser’ (default), ‘notebook’ and ‘png’.

  • **kwargs (dict[str, object]) – Additional keyword arguments to pass to the scatter_3d function of the plotly.express package.

Raises:

ValueError – If the cage is not loaded.

Warns:
  • UserWarning – If the cavity is not detected, a warning is issued.

  • UserWarning – If the cluster is not packed, a warning is issued.

  • UserWarning – If the openings are not detected, a warning is issued.

class AtomPacker.Cavity(grid: ndarray, step: float, probe_in: float, probe_out: float, removal_distance: float, volume_cutoff: float, vertices: ndarray, surface: str = 'SES')[source]

A class representing a cavity in a molecular structure.

The Cavity class provides methods for detecting cavities, calculating cavity properties, and exporting cavity data.

property coordinates: ndarray

Get the xyz coordinates of the cavity.

detect_openings(openings_cutoff: int = 1, verbose: bool = False) None[source]

Detect openings in the cavity.

This method detects openings in the cavity represented by the Cavity object. It uses a grid-based algorithm to identify openings based on the depth of the cavity points.

Parameters:
  • openings_cutoff (float, optional) – The cutoff value for detecting openings (default is 1). The minimum number of points in an opening to be considered valid.

  • verbose (bool, optional) – If True, print detailed information during processing (default is False).

Returns:

An AtomPacker.Openings object containing the detected openings.

Return type:

Openings

preview(renderer: str = 'browser', opacity: float = 1.0, **kwargs: dict[str, object]) None[source]

Preview the cavity in a 3D viewer.

Parameters:
  • renderer (str, optional) – The renderer to use for the 3D viewer. Supported renderers are ‘browser’ (default), ‘notebook’ and ‘png’.

  • opacity (float, optional) – The opacity of the atoms in the 3D viewer, by default 1.0. The opacity value ranges from 0.0 (completely transparent) to 1.0 (completely opaque).

  • **kwargs (dict[str, object], optional) – Additional keyword arguments to pass to the scatter_3d function of the plotly.express package.

save(filename: str = 'cavity.pdb') None[source]

Export the cavity data to a file.

This method exports the cavity data represented by the Cavity object to a file in a specified format. The supported formats include PDB, and XYZ.

Parameters:

filename (str) – The name of the file to export the cavity data to. The filename of the structure file. The file format is determined by the suffix. Supported formats are: .pdb, .xyz.

select_cavity(indices: list[int]) None[source]

Select a cavity from the list of detected cavities.

This method selects a cavity from the list of detected cavities based on the indices of the cavity in the list.

Parameters:

indices (list[int]) –

The indices of the cavities to select. Cavity points in the 3D grid (grid[nx][ny][nz]). Cavities array has integer labels in each position, that are:

  • -1: solvent points;

  • 0: cage points;

  • 1: empty space points;

  • >=2: cavity points.

property volume: float64

Calculate the volume of the cavity.

This method calculates the volume of the cavity represented by the Cavity object. It uses a grid-based algorithm to estimate the volume of the cavity based on the number of grid points inside the cavity.

class AtomPacker.Cluster(cluster: Cluster, cavity: Cavity, log: DataFrame | None = None)[source]

A class representing a cluster of atoms inside a supramolecular cage.

The Cluster class provides methods for calculating cluster properties and exporting cluster data.

property coordinates: ndarray

Return the coordinates of the atoms in the cluster.

diameter(method: str = 'maximum') float[source]

Calculate the diameter of the cluster.

Parameters:

method (str, optional) –

The method to use for calculating the diameter of the cluster. Supported methods are ‘maximum’ (default), ‘shape’, and ‘volume’.

  • ’maximum’: calculates the maximum distance between any two atoms in the cluster, considering their radius.

  • ’shape’: uses ase.cluster.Cluster.get_diameter(‘shape’) to return the averaged diameter calculated from the directions given by the defined surfaces.

  • ’volume’: uses ase.cluster.Cluster.get_diameter(‘volume’) to return the diameter of a sphere with the same volume as the atoms.

Returns:

The diameter of the cluster.

Return type:

float

property maximum_number_of_atoms: int

Calculate the maximum number of atoms in the cluster.

The maximum number of atoms is calculated based on the volume of the cavity divided by the volume of the atom in the cluster. The formula is given by:

\[N_{max} = \left \lceil \frac{V_{cav}}{V_a} \right \rceil = \left \lceil \frac{V_{cav}}{\frac{4}{3}\pi R_a^3} \ \right \rceil\]

where \(V_{cav}\) is the volume of the cavity, \(V_a\) is the volume of the atom in the cluster, \(R_a\) is the radius of the atom in the cluster, and \(\lceil \cdot \rceil\) is the ceiling function.

property number_of_atoms: int

Get the number of atoms in the cluster.

preview(renderer: str = 'browser', opacity: float = 1.0, **kwargs: dict[str, object]) None[source]

Preview the cavity in a 3D viewer.

Parameters:
  • renderer (str, optional) – The renderer to use for the 3D viewer. Supported renderers are ‘browser’ (default), ‘notebook’ and ‘png’.

  • opacity (float, optional) – The opacity of the atoms in the 3D viewer, by default 1.0. The opacity value ranges from 0.0 (completely transparent) to 1.0 (completely opaque).

  • **kwargs (dict[str, object]) – Additional keyword arguments to pass to the scatter_3d function of the plotly.express package.

save(filename: str = 'cluster.pdb') None[source]

Export the cluster data to a file.

This method exports the cluster data represented by the Cluster object to a file in a specified format. The supported formats include PDB, and XYZ.

Parameters:

filename (str, optional) – The name of the file to export the cavity data to, by default “cluster.pdb”. The file format is determined by the suffix. Supported formats are: .pdb, .xyz.

property summary: DataFrame

Print a summary of the cluster properties.

class AtomPacker.Openings(cavities: ndarray, step: float, vertices: ndarray, openings_cutoff: float = 1, verbose: bool = False)[source]

A class representing openings in a cavity structure.

The Openings class provides methods for detecting openings, calculating their properties, and exporting openings data.

property coordinates: ndarray

Get the xyz coordinates of the cavity.

preview(renderer: str = 'browser', opacity: float = 1.0, **kwargs: dict[str, object]) None[source]

Preview the openings in a 3D viewer.

Parameters:
  • renderer (str, optional) – The renderer to use for the 3D viewer. Supported renderers are ‘browser’ (default), ‘notebook’ and ‘png’.

  • opacity (float, optional) – The opacity of the atoms in the 3D viewer, by default 1.0. The opacity value ranges from 0.0 (completely transparent) to 1.0 (completely opaque).

  • **kwargs (dict[str, object], optional) – Additional keyword arguments to pass to the scatter_3d function of the plotly.express package.

save(filename: str = 'openings.pdb') None[source]

Save the openings to a PDB file.

Parameters:

filename (str, optional) – The name of the output PDB file (default is “openings.pdb”).

AtomPacker.get_coordinates(grid: ndarray, step: float, vertices: ndarray) ndarray[source]

Convert a grid of points to a numpy.ndarray of xyz coordinates.

Parameters:
  • grid (numpy.ndarray) – The grid points.

  • step (float) – The step size of the grid.

  • vertices (numpy.ndarray) – The vertices (origin, X-axis, Y-axis, Z-axis) of the grid.

Returns:

coordinates – The xyz coordinates of the grid points.

Return type:

numpy.ndarray

AtomPacker.get_depths(grid: ndarray, step: float) ndarray[source]

Get the depth of a grid of cavity points.

Parameters:
  • grid (numpy.ndarray) – The grid of cavity points.

  • step (float) – The step size of the grid.

Returns:

depths – The depths of the cavity points.

Return type:

numpy.ndarray

AtomPacker.get_lattice_constants(atom_type: str, lattice_type: str) tuple[float, float, float] | tuple[float, float] | float | None[source]

Get the lattice constants for a given atom type and lattice type.

Parameters:
  • atom_type (str) – The atomic symbol of the atom.

  • lattice_type (str) – The type of lattice in the cluster. The available lattice types are ‘bcc’, ‘fcc’, ‘hcp’ and ‘sc’.

Returns:

The lattice constants for the given atom type and lattice type. If None, use ase.data.reference_states.

Return type:

tuple[float, float, float] | tuple[float, float] | float | None

Note

The lattice constants will be fetched from AtomPacker.data.lattice_constants if available. If not, the experimental values from ase.data will be used.

Raises:

ValueError – If the lattice type is not valid.

AtomPacker.load_mmcif(filename: str, vdw: Dict[str, float] | None = None) Universe[source]

Load a PDBx/mmCIF file into the MDAnalysis.Univese object.

Parameters:
  • filename (str) – The filename of the structure file.

  • vdw (Dict[str, float]) – A dictionary containing the van der Waals radii for each atom type, by default None. If None, the radii will be looked up from the pyKVFinder package.

Returns:

universe – The MDAnalysis universe object containing the structure.

Return type:

MDAnalysis.Univese

AtomPacker.load_mol2(filename: str, vdw: Dict[str, float] | None = None) Universe[source]

Load a MOL2 file into the MDAnalysis.Univese object.

Parameters:
  • filename (str) – The filename of the structure file.

  • vdw (Dict[str, float]) – A dictionary containing the van der Waals radii for each atom type, by default None. If None, the radii will be looked up from the pyKVFinder package.

Returns:

universe – The MDAnalysis universe object containing the structure.

Return type:

MDAnalysis.Univese

AtomPacker.load_pdb(filename: str, vdw: Dict[str, float] | None = None) Universe[source]

Load a PDB file into the MDAnalysis.Univese object.

Parameters:
  • filename (str) – The filename of the structure file.

  • vdw (Dict[str, float]) – A dictionary containing the van der Waals radii for each atom type, by default None. If None, the radii will be looked up from the pyKVFinder package.

Returns:

universe – The MDAnalysis universe object containing the structure.

Return type:

MDAnalysis.Univese

AtomPacker.load_xyz(filename: str, vdw: Dict[str, float] | None = None) Universe[source]

Load a XYZ file into the MDAnalysis.Univese object.

Parameters:
  • filename (str) – The filename of the structure file.

  • vdw (Dict[str, float]) – A dictionary containing the van der Waals radii for each atom type, by default None. If None, the radii will be looked up from the pyKVFinder package.

Returns:

universe – The MDAnalysis universe object containing the structure.

Return type:

MDAnalysis.Univese