πŸ“– API Documentation

moxel package

MOXΡλ is a Python package for parallel calculation of energy voxels, with emphasis on reticular chemistry.

Note

Currently, interactions are modelled with the Lennard-Jones (LJ) potential.

Submodules

moxel.utils module

class moxel.utils.Grid(grid_size=25, cutoff=10, epsilon=50, sigma=2.5)

Bases: object

A 3D energy grid over a crystal structure.

Parameters:
  • grid_size (int, default=25) – Number of grid points along each dimension.

  • cutoff (float, default=10) – Cutoff radius (β„«) for the LJ potential.

  • epsilon (float, default=50) – Epsilon value (Ξ΅/K) of the probe atom.

  • sigma (float, default=2.5) – Sigma value (Οƒ/β„«) of the probe atom.

structure

Available only after Grid.load_structure() has been called.

Type:

pymatgen.core.structure.Structure

structure_name

Available only after Grid.load_structure() has been called.

Type:

str

cubic_box

Available only after Grid.calculate() has been called.

Type:

bool

voxels

Available only after Grid.calculate() has been called.

Type:

array of shape (grid_size, grid_size, grid_size)

calculate(cubic_box=False, length=30, potential='lj')

Iterate over the grid and return voxels.

For computational efficiency and to assure (approximately) the same spatial resolution, the grid is overlayed over a supercell scaled according to MIC, see mic_scale_factors().

If lattice angles are significantly different than 90Β°, to avoid distortions set cubic_box to True. In this case, the grid is overlayed over a cubic box of size length centered at the origin but periodicity is no longer guaranteed.

Parameters:
  • potential (str, default='lj') – The potential used to calculate voxels. Currently, only the LJ potential is supported.

  • cubic_box (bool, default=False) – If True, the simulation box is cubic.

  • length (float, default=30) – The size of the cubic box in Γ…. Takes effect only if cubic_box=True.

Returns:

voxels – The energy voxels as \(e^{-\beta \mathcal{V}}\), to ensure numerical stability.

Return type:

array of shape (grid_size, grid_size, grid_size)

Notes

For structures that can not be processsed, their voxels are filled with zeros.

lj_potential(coords)

Calculate LJ potential at cartesian or fractional coordinates.

Parameters:

coordinates (array_like of shape (3,)) – If cubic_box=True cartesian. Else, fractional.

Returns:

energy – Energy as \(e^{-\beta \mathcal{V}}\), to ensure numerical stability.

Return type:

float

load_structure(cif_pathname)

Load a crystal structure from a .cif file.

Parameters:

cif_pathname (str) – Pathname to the .cif file.

moxel.utils.batch_clean_and_merge(batch_dirnames, out_pathname=None)

Clean a single batch, or first clean and then merge multiple batches. All batches must have the form:

batch
β”œβ”€β”€voxels.npy
└──names.json

Cleaning is required since the voxels for some structures might be zero, see Grid.calculate().

If len(batch_dirnames) == 1 the cleaned voxels for are stored under batch_dirnames[0]/clean_voxels.npy and the names of their corresponding structures under batch_dirnames[0]/clean_names.json.

If len(batch_dirnames) > 1 the voxels (cleaned and merged) are stored under out_pathname/clean_voxels.npy and the names of their corresponding structures under out_pathname/clean_names.json. That is:

out_pathname
β”œβ”€β”€clean_voxels.npy
└──clean_names.json
Parameters:
  • batch_dirnames (list) – List of pathnames to the directories of the batches.

  • out_pathname (str, optional) – Pathname to the directory holding the clean voxels and names. The directory is created if it doesn’t exist. Takes effect only if len(batch_dirnames) > 1. If None voxels and names are stored under ./clean_voxels.npy and ./clean_names.json.

Returns:

exit_status – If no voxels are missing 0 else 1.

Return type:

int

moxel.utils.mic_scale_factors(r, lattice_vectors)

Return scale factors to satisfy minimum image convention (MIC) [1] .

Parameters:
  • r (float) – The cutoff radius used in MIC convetion.

  • lattice_vectors (array of shape (3, 3)) – The lattice vectors of the unit cell. Each row corresponds to a lattice vector.

Returns:

scale_factors – scale_factors[i] scales lattice_vectors[i].

Return type:

array of shape (3,)

References

moxel.utils.plot_voxels(voxels, *, fill_pattern=None, colorbar=True, cmap='viridis', **kwargs)

Visualizing voxels with Axes3d.voxels.

Parameters:
  • voxels (3D array) –

  • fill_pattern (3D array of bool, optional) – A 3D array of truthy values, indicating which voxels to fill. If not specified, all voxels are filled.

  • colorbar (bool, default=True) – Whether to include a colorbar.

  • cmap (str, default='viridis') – Colormap that colorizes the voxels based on their value.

  • **kwargs –

    Valid keyword arguments for Axes3d.voxels.

    Warning

    Do not pass the argument facecolors of Axes3d.voxels. This argument is used under the hood by plot_voxels_mpl() to generate the colors of the voxels based on the specified cmap.

Returns:

fig

Return type:

Figure

moxel.utils.voxels_from_dir(cif_dirname, grid_size=25, cutoff=10, epsilon=50, sigma=2.5, cubic_box=False, length=30, out_pathname=None)

Calculate voxels from .cif files in cif_dirname and save them in out_pathname as array of shape (n_samples, grid_size, grid_size, grid_size), where n_samples is the number of .cif files in cif_dirname.

Parameters:
  • cif_dirname (str) – Pathname to the directory containing the .cif files.

  • grid_size (int, default=25) – Number of grid points along each dimension.

  • cutoff (float, default=10) – Cutoff radius (β„«) for the LJ potential.

  • epsilon (float, default=50) – Epsilon value (Ξ΅/K) of the probe atom.

  • cubic_box (bool, default=False) – If True, the simulation box is cubic.

  • length (float, default=30) – The size of the cubic box in Γ…. Takes effect only if cubic_box=True.

  • sigma (float, default=25) – Sigma value (Οƒ/β„«) of the probe atom.

  • out_pathname (str, optional) – Pathname to the file holding the voxels. If not specified, voxels are stored in ./voxels.npy.

Notes

  • Samples in output array follow the order in sorted(os.listdir(cif_dirname)).

  • For structures that can not be processsed, their voxels are filled with zeros.

moxel.utils.voxels_from_file(cif_pathname, grid_size=25, cutoff=10, epsilon=50, sigma=2.5, cubic_box=False, length=30, only_voxels=True)

Return voxels from .cif file.

Parameters:
  • cif_pathname (str) – Pathname to the .cif file.

  • grid_size (int, default=25) – Number of grid points along each dimension.

  • cutoff (float, default=10) – Cutoff radius (β„«) for the LJ potential.

  • epsilon (float, default=50) – Epsilon value (Ξ΅/K) of the probe atom.

  • sigma (float, default=25) – Sigma value (Οƒ/β„«) of the probe atom.

  • cubic_box (bool, default=False) – If True, the simulation box is cubic.

  • length (float, default=30) – The size of the cubic box in Γ…. Takes effect only if cubic_box=True.

  • only_voxels (bool, default=True) – Determines out type.

Returns:

out – If only_voxels=True, array of shape (grid_size, grid_size, grid_size). Otherwise, Grid.

Return type:

array or Grid

Notes

  • For structures that can not be processsed, their voxels are filled with zeros.

moxel.utils.voxels_from_files(cif_pathnames, grid_size=25, cutoff=10, epsilon=50, sigma=2.5, cubic_box=False, length=30, out_pathname=None)

Calculate voxels from a list of .cif files and save them in out_pathname as array of shape (n_samples, grid_size, grid_size, grid_size), where n_samples is the number of is the number of .cif files in cif_pathnames.

Parameters:
  • cif_pathanmes (list) – List of pathnames to the .cif files.

  • grid_size (int, default=25) – Number of grid points along each dimension.

  • cutoff (float, default=10) – Cutoff radius (β„«) for the LJ potential.

  • epsilon (float, default=50) – Epsilon value (Ξ΅/K) of the probe atom.

  • sigma (float, default=25) – Sigma value (Οƒ/β„«) of the probe atom.

  • cubic_box (bool, default=False) – If True, the simulation box is cubic.

  • length (float, default=30) – The size of the cubic box in Γ…. Takes effect only if cubic_box=True.

  • out_pathname (str, optional) – Pathname to the file holding the voxels. If not specified, voxels are stored in ./voxels.npy.

Notes

  • Samples in output array follow the order in cif_pathnames.

  • For structures that can not be processsed, their voxels are filled with zeros.