π API Documentationο
moxel packageο
MOXΡλ is a Python package for parallel calculation of energy voxels, with emphasis on reticular chemistry.
Submodulesο
moxel.utilsο
This module provides helper functions for creating voxels.
Note
Currently, interactions are modelled with the Lennard-Jones (LJ) potential.
- class moxel.utils.Grid(grid_size=25, cutoff=10, epsilon=50, sigma=2.5)[source]ο
Bases:
objectA 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,)*3
- calculate(cubic_box=False, length=30, potential='lj', n_jobs=None)[source]ο
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_boxtoTrue. In this case, the grid is overlayed over a cubic box of sizelengthcentered 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.n_jobs (int, optional) β Number of jobs to run in parallel. If
None, all processors are used.
- Returns:
voxels β The energy voxels as \(e^{-\beta \mathcal{V}}\), to ensure numerical stability.
- Return type:
array of shape (grid_size,)*3
Notes
For structures that can not be processsed, their voxels are filled with zeros.
- lj_potential(coords)[source]ο
Calculate LJ potential at cartesian or fractional coordinates.
- Parameters:
coordinates (array_like of shape (3,)) β If
cubic_box == Truecartesian. Else, fractional.- Returns:
energy β Energy as \(e^{-\beta \mathcal{V}}\), to ensure numerical stability.
- Return type:
float
- moxel.utils.batch_clean(batch_dirname)[source]ο
Clean a single batch.
The batch must have the form:
batch βββvoxels.npy βββnames.json
Cleaning is required since the voxels for some structures might be zero, see
Grid.calculate(). After cleaning, the directory has the form:batch βββvoxels.npy βββnames.json βββclean_voxels.npy βββclean_names.json
- Parameters:
batch_dirname (str) β Pathname to the directory which requires cleaning.
- Returns:
exit_status β If no voxels are missing
0else1.- Return type:
int
- moxel.utils.get_names(fname)[source]ο
Load a list of material names saved in
.jsonformat.- Parameters:
fname (str) β Pathname to the
.jsonfile.- Returns:
names
- Return type:
list
- moxel.utils.mic_scale_factors(r, lattice_vectors)[source]ο
Return scale factors to satisfy minimum image convention [MIC].
- 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]scaleslattice_vectors[i].- Return type:
array of shape (3,)
References
[MIC]Smith, βThe Minimum Image Convention in Non-Cubic MD Cellsβ, 1989.
- moxel.utils.voxels_from_dir(cif_dirname, out_pathname, grid_size=25, cutoff=10, epsilon=50, sigma=2.5, cubic_box=False, length=30, n_jobs=None)[source]ο
Calculate voxels from a directory of
.ciffiles and save them underout_pathnameasnumpy.arrayof shape(n_samples, grid_size, grid_size, grid_size), wheren_samples == len(cif_pathnames).After processing the following files are created:
out_pathname βββvoxels.npy βββnames.jsonThe file
names.jsonstores the names of the materials, which might be useful for later indexing.- Parameters:
cif_dirname (str) β Pathname to the directory containing the
.ciffiles.out_pathname (str) β Pathname to the directory under which voxels are stored.
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.n_jobs (int, optional) β Number of jobs to run in parallel. If
None, all processors are used.
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, n_jobs=None, only_voxels=True)[source]ο
Return voxels from
.ciffile.- Parameters:
cif_pathname (str) β Pathname to the
.ciffile.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.n_jobs (int, optional) β Number of jobs to run in parallel. If
None, all processors are used.only_voxels (bool, default=True) β Determines
outtype.
- Returns:
out β If
only_voxels == True, array of shape(grid_size,)*3. Otherwise,Grid.- Return type:
arrayorGrid
Notes
For structures that can not be processsed, their voxels are filled with zeros.
- moxel.utils.voxels_from_files(cif_pathnames, out_pathname, grid_size=25, cutoff=10, epsilon=50, sigma=2.5, cubic_box=False, length=30, n_jobs=None)[source]ο
Calculate voxels from a list of
.ciffiles and store them underout_pathnameasnumpy.arrayof shape(n_samples, grid_size, grid_size, grid_size), wheren_samples == len(cif_pathnames).After processing the following files are created:
out_pathname βββvoxels.npy βββnames.jsonThe file
names.jsonstores the names of the materials, which might be useful for later indexing.- Parameters:
cif_pathnames (list) β List of pathnames to the
.ciffiles.out_pathname (str) β Pathname to the directory under which voxels are stored.
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.n_jobs (int, optional) β Number of jobs to run in parallel. If
None, all processors are used.
Notes
Samples in output array follow the order in
cif_pathnames.For structures that can not be processsed, their voxels are filled with zeros.
moxel.visualizeο
This module provides helper functions for visualizing voxels.
Tip
For faster rendering you should prefer plot_voxels_pv().
- moxel.visualize.plot_voxels_mpl(voxels, *, fill_pattern=None, colorbar=True, cmap='viridis', **kwargs)[source]ο
Visualize 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.
Note
Do not pass the argument
facecolorsof Axes3d.voxels. This argument is used under the hood byplot_voxels()to generate the colors of the voxels based on the specifiedcmap.
- Returns:
fig
- Return type:
Examples
>>> fig = plot_voxels_mpl(voxels, cmap='coolwarm')
- moxel.visualize.plot_voxels_pv(voxels, **kwargs)[source]ο
Visualize voxels with Plotter.add_volume
Note
For interactive plots in Jupyter:
pip install "pyvista[jupyter]"- Parameters:
voxels (3D array)
**kwargs (dict, optional) β Valid keyword arguments for Plotter.add_volume.
Examples
>>> plot_voxels_pv(voxels, cmap='coolwarm')