pywindow package

Submodules

pywindow.io_tools module

Module contains classes for input/output processing.

class pywindow.io_tools.Input[source]

Bases: object

Class used to load and process input files.

load_file(filepath)[source]

This function opens any type of a readable file and decompose the file object into a list, for each line, of lists containing splitted line strings using space as a spacer.

Parameters:filepath (str) – The full path or a relative path to any type of file.
Returns:Returns a dictionary containing the molecular information extracted from the input files. This information will vary with file type and information stored in it. The data is sorted into lists that contain one feature for example key atom_id: [atom_id_1, atom_id_2] Over the process of analysis this dictionary will be updated with new data.
Return type:dict
load_rdkit_mol(mol)[source]

Return molecular data from rdkit.Chem.rdchem.Mol object.

Parameters:mol (rdkit.Chem.rdchem.Mol) – A molecule object from RDKit.
Returns:A dictionary with elements and coordinates as keys containing molecular data extracted from rdkit.Chem.rdchem.Mol object.
Return type:dict
class pywindow.io_tools.Output[source]

Bases: object

Class used to process and save output files.

dump2file(obj, filepath, override=False, **kwargs)[source]

Dump a dictionary into a file. (Extensions: XYZ or PDB)

Parameters:
  • obj (dict) – A dictionary containing molecular information.
  • filepath (str) – The filepath for the dumped file.
  • override (bool) – If True, any file in the filepath will be override. (default=False)
dump2json(obj, filepath, override=False, **kwargs)[source]

Dump a dictionary into a JSON dictionary.

Uses the json.dump() function.

Parameters:
  • obj (dict) – A dictionary to be dumpped as JSON file.
  • filepath (str) – The filepath for the dumped file.
  • override (bool) – If True, any file in the filepath will be override. (default=False)

pywindow.molecular module

Defines MolecularSystem and Molecule classes.

This module is the most important part of the pywindow package, as it is at the frontfront of the interaction with the user. The two main classes defined here: MolecularSystem and Molecule are used to store and analyse single molecules or assemblies of single molecules.

The MolecularSystem is used as a first step to the analysis. It allows to load data, to refine it (rebuild molecules in a periodic system, decipher force field atom ids) and to extract single molecules for analysis as Molecule instances.

To get started see MolecularSystem.

To get started with the analysis of Molecular Dynamic trajectories go to pywindow.trajectory.

class pywindow.molecular.MolecularSystem[source]

Bases: object

Container for the molecular system.

To load input and initialise MolecularSystem, one of the MolecularSystem classmethods (load_file(), load_rdkit_mol() or load_system()) should be used. MolecularSystem should not be initialised by itself.

Examples

  1. Using file as an input:
pywindow.molecular.MolecularSystem.load_file(`filepath`)
  1. Using RDKit molecule object as an input:
pywindow.molecular.MolecularSystem.load_rdkit_mol(rdkit.Chem.rdchem.Mol)
  1. Using a dictionary (or another MoleculeSystem.system) as input:
pywindow.molecular.MolecularSystem.load_system({...})
system_id

str or int – The input filename or user defined.

system

dict – A dictionary containing all the information extracted from input.

molecules

list – A list containing all the returned Molecule s after using make_modular().

decipher_atom_keys(forcefield='DLF', dict_key='atom_ids')[source]

Decipher force field atom ids.

This takes all values in MolecularSystem.system['atom_ids'] that match force field type criteria and creates MolecularSystem.system['elements'] with the corresponding periodic table of elements equivalents.

If a forcefield is not supported by this method, the MolecularSystem.swap_atom_keys() can be used instead.

DLF stands for DL_F notation.

See: C. W. Yong, Descriptions and Implementations of DL_F Notation: A Natural Chemical Expression System of Atom Types for Molecular Simulations, J. Chem. Inf. Model., 2016, 56, 1405–1409.

Parameters:
  • forcefield (str) – The forcefield used to decipher atom ids. Allowed (not case sensitive): ‘OPLS’, ‘OPLS2005’, ‘OPLSAA’, ‘OPLS3’, ‘DLF’, ‘DL_F’. (default=’DLF’)
  • dict_key (str) – The MolecularSystem.system dictionary key to the array containing the force field atom ids. (default=’atom_ids’)
Returns:

None

Return type:

NoneType

dump_system(filepath=None, modular=False, **kwargs)[source]

Dump a MolecularSystem to a file (PDB or XYZ).

Kwargs are passed to pywindow.io_tools.Output.dump2file().

Parameters:
Returns:

None

Return type:

NoneType

dump_system_json(filepath=None, modular=False, **kwargs)[source]

Dump a MolecularSystem to a JSON dictionary.

The dumped JSON dictionary, with MolecularSystem, can then be loaded through a JSON loader and then through load_system() to retrieve a MolecularSystem.

Kwargs are passed to pywindow.io_tools.Output.dump2json().

Parameters:
Returns:

None

Return type:

NoneType

classmethod load_file(filepath)[source]

Create a MolecularSystem from an input file.

Recognized input file formats: XYZ, PDB and MOL (V3000).

Parameters:filepath (str) – The input’s filepath.
Returns:MolecularSystem
Return type:pywindow.molecular.MolecularSystem
classmethod load_rdkit_mol(mol)[source]

Create a MolecularSystem from rdkit.Chem.rdchem.Mol.

Parameters:mol (rdkit.Chem.rdchem.Mol) – An RDKit molecule object.
Returns:MolecularSystem
Return type:pywindow.molecular.MolecularSystem
classmethod load_system(dict_, system_id='system')[source]

Create a MolecularSystem from a python dict.

As the loaded MolecularSystem is storred as a dict in the MolecularSystem.system it can also be loaded directly from a dict input. This feature is used by trajectory that extracts trajectory frames as dictionaries and returns them as MolecularSystem objects through this classmethod.

Parameters:
  • dict (dict) – A python dictionary.
  • system_id (str or :class:’int’, optional) – Inherited or user defined system id. (default=’system’)
Returns:

MolecularSystem

Return type:

pywindow.molecular.MolecularSystem

make_modular(rebuild=False)[source]

Find and return all Molecule s in MolecularSystem.

This function populates MolecularSystem.molecules with Molecule s.

Parameters:rebuild (bool) – If True, run first the rebuild_system(). (default=False)
Returns:None
Return type:NoneType
rebuild_system(override=False, **kwargs)[source]

Rebuild molecules in molecular system.

Parameters:override (bool, optional (default=False)) – If False the rebuild molecular system is returned as a new MolecularSystem, if True, the current MolecularSystem is modified.
swap_atom_keys(swap_dict, dict_key='atom_ids')[source]

Swap a force field atom id for another user-defined value.

This modified all values in MolecularSystem.system['atom_ids'] that match criteria.

This function can be used to decipher a whole forcefield if an appropriate dictionary is passed to the function.

Example

In this example all atom ids ‘he’ will be exchanged to ‘H’.

pywindow.molecular.MolecularSystem.swap_atom_keys({'he': 'H'})
Parameters:
  • swap_dict (dict) – A dictionary containg force field atom ids (keys) to be swapped with corresponding values (keys’ arguments).
  • dict_key (str) – A key in MolecularSystem.system dictionary to perform the atom keys swapping operation on. (default=’atom_ids’)
Returns:

None

Return type:

NoneType

system_to_molecule()[source]

Return MolecularSystem as a Molecule directly.

Only to be used conditionally, when the MolecularSystem is a discrete molecule and no input pre-processing is required.

Returns:Molecule
Return type:pywindow.molecular.Molecule
class pywindow.molecular.Molecule(mol, system_name, mol_id)[source]

Bases: pywindow.molecular._Shape

Container for a single molecule.

This class is meant for the analysis of single molecules, molecular pores especially. The object passed to this class should therefore be a finite and interconnected individuum.

This class should not be initialised directly, but result from MolecularSystem.system_to_molecule() or MolecularSystem.make_modular().

Methods in Molecule allow to calculate:

  1. The maximum diameter of a molecule.
  2. The average diameter of a molecule.
  3. The intrinsic void diameter of a molecule.
  4. The intrinsic void volume of a molecule.
  5. The optimised intrinsic void diameter of a molecule.
  6. The optimised intrinsic void volume of a molecule.
  7. The circular diameter of a window of a molecule.
mol

dict – The Molecular.System.system dictionary passed to the Molecule which is esentially a container of the information that compose a molecular entity, such as the coordinates and atom ids and/or elements.

no_of_atoms

int – The number of atoms in the molecule.

elements

numpy.array – An array containing the elements, as strings, composing the molecule.

atom_ids

numpy.array (conditional) – If the Molecule.mol contains ‘atom_ids’ keyword, the force field ids of the elements.

coordinates

numpy.array – The x, y and z atomic Cartesian coordinates of all elements.

parent_system

str – The name of MolecularSystem passed to Molecule.

molecule_id

any – The molecule id passed when initialising Molecule.

properties

dict – A dictionary that is populated by the output of Molecule methods.

calculate_average_diameter(**kwargs)[source]

Return the average diamension of a molecule.

Returns:The average dimension of the molecule.
Return type:float
calculate_centre_of_mass()[source]

Return the xyz coordinates of the centre of mass of a molecule.

Returns:The centre of mass of the molecule.
Return type:numpy.array
calculate_maximum_diameter()[source]

Return the maximum diamension of a molecule.

Returns:The maximum dimension of the molecule.
Return type:float
calculate_pore_diameter()[source]

Return the intrinsic pore diameter.

Returns:The intrinsic pore diameter.
Return type:float
calculate_pore_diameter_opt(**kwargs)[source]

Return the intrinsic pore diameter (for the optimised pore centre).

Similarly to calculate_pore_diameter() this method returns the the intrinsic pore diameter, however, first a better approximation of the pore centre is found with optimisation.

Returns:The intrinsic pore diameter.
Return type:float
calculate_pore_volume()[source]

Return the intrinsic pore volume.

Returns:The intrinsic pore volume.
Return type:float
calculate_pore_volume_opt(**kwargs)[source]

Return the intrinsic pore volume (for the optimised pore centre).

Similarly to calculate_pore_volume() this method returns the the volume intrinsic pore diameter, however, for the calculate_pore_diameter_opt() returned value.

Returns:The intrinsic pore volume.
Return type:float
calculate_windows(**kwargs)[source]

Return the diameters of all windows in a molecule.

This function first finds and then measures the diameters of all the window in the molecule.

Returns:
  • numpy.array – An array of windows’ diameters.
  • NoneType – If no windows were found.
dump_molecule(filepath=None, include_coms=False, **kwargs)[source]

Dump a Molecule to a file (PDB or XYZ).

Kwargs are passed to pywindow.io_tools.Output.dump2file().

For validation purposes an overlay of window centres and COMs can also be dumped as:

He - for the centre of mass

Ne - for the centre of the optimised cavity

Ar - for the centres of each found window

Parameters:
  • filepath (str) – The filepath for the dumped file. If None, the file is dumped localy with molecule_id as filename. (defualt=None)
  • include_coms (bool) – If True, dump also with an overlay of window centres and COMs. (default=False)
Returns:

None

Return type:

NoneType

dump_properties_json(filepath=None, molecular=False, **kwargs)[source]

Dump content of Molecule.properties to a JSON dictionary.

Parameters:
  • filepath (str) – The filepath for the dumped file. If None, the file is dumped localy with molecule_id as filename. (defualt=None)
  • molecular (bool) – If False, dump only the content of Molecule.properties, if True, dump all the information about Molecule.
Returns:

None

Return type:

NoneType

full_analysis(ncpus=1, **kwargs)[source]

Perform a full structural analysis of a molecule.

This invokes other methods:

Parameters:ncpus (int) – Number of CPUs used for the parallelised parts of pywindow.utilities.find_windows(). (default=1=serial)
Returns:The updated Molecule.properties with returns of all used methods.
Return type:Molecule.properties
molecular_weight()[source]

Return the molecular weight of a molecule.

Returns:The molecular weight of the molecule.
Return type:float
shift_to_origin(**kwargs)[source]

Shift a molecule to Origin.

This function takes the molecule’s coordinates and adjust them so that the centre of mass of the molecule coincides with the origin of the coordinate system.

Returns:None
Return type:NoneType

pywindow.tables module

This module containes general-purpose chemical data (aka tables).

Sources:

1. www.ccdc.cam.ac.uk/Lists/ResourceFileList/Elemental_Radii.xlsx, (access date: 13 Oct 2015)

2. C. W. Yong, ‘DL_FIELD - A force field and model development tool for DL_POLY’, R. Blake, Ed., CSE Frontier, STFC Computational Science and Engineering, Daresbury Laboratory, UK, p38-40 (2010)

3. https://www.ccdc.cam.ac.uk/support-and-resources/ccdcresources/ Elemental_Radii.xlsx (access date: 26 Jul 2017)

Note:

Added atomic mass, atomic van der Waals radius and covalent radius, that equals 1, for a dummy atom X. This is for the shape descriptors calculations to use with functions for molecular shape descriptors that require mass.

pywindow.trajectory module

Module intended for the analysis of molecular dynamics trajectories.

The trajectory file (DL_POLY_C:HISTORY, PDB or XYZ) should be loaded with the one of the corresponding classes (DLPOLY, PDB or XYZ, respectively).

Example

In this example a DL_POLY_C HISTORY trajectory file is loaded.

pywindow.trajectory.DLPOLY('path/to/HISTORY')

Then, each of the trajectory frames can be extracted and returned as a pywindow.molecular.MolecularSystem object for analysis. See pywindow.molecular docstring for more information.

Alternatively, the analysis can be performed on a whole or a chunk of the trajectory with the analysis() function. The benefit is that the analysis can be performed in parallel and the results stored as a single JSON dictionary in a straightforward way. Also, the deciphering of the force field atom ids and the rebuilding of molecules can be applied to each frame in a consitent and automated manner. The downfall is that at the moment it is not possible to choose the set of parameters that are being calculated in the pywindow.molecular.Molecule as the pywindow.molecular.Molecule.full_analysis() is invoked by default. However, the computational cost of calculating majority of the structural properties is miniscule and it is usually the pywindow.molecular.MolecularSystem.rebuild_system() step that is the bottleneck.

class pywindow.trajectory.DLPOLY(filepath)[source]

Bases: object

A container for a DL_POLY_C type trajectory (HISTORY).

This function takes a DL_POLY_C trajectory file and maps it for the binary points in the file where each frame starts/ends. This way the process is fast, as it not require loading the trajectory into computer memory. When a frame is being extracted, it is only this frame that gets loaded to the memory.

Frames can be accessed individually and loaded as an unmodified string, returned as a pywindow.molecular.MolecularSystem (and analysed), dumped as PDB or XYZ or JSON (if dumped as a pywindow.molecular.MolecularSystem.system)

filepath

str – The filepath.

system_id

str – The system id inherited from the filename.

frames

dict – A dictionary that is populated, on the fly, with the extracted frames.

analysis_output

dict – A dictionary that is populated, on the fly, with the analysis output.

analysis(frames='all', ncpus=1, _ncpus=1, override=False, **kwargs)[source]

Perform structural analysis on a frame/ set of frames.

Depending on the passed parameters a frame, a list of particular frames, a range of frames (from, to), or all frames can be analysed with this function.

The analysis is performed on each frame and each discrete molecule in that frame separately. The steps are as follows:

  1. A frame is extracted and returned as a MolecularSystem.
  2. If swap_atoms is set the atom ids are swapped.
  3. If forcefield is set the atom ids are deciphered.
  4. If rebuild is set the molecules in the system are rebuild.
  5. Each discrete molecule is extracted as Molecule
  6. Each molecule is analysed with Molecule.full_analysis()
  7. Analysis output populates the analysis_output dictionary.

As the analysis of trajectories often have to be unique, many options are conditional.

A side effect of this function is that the analysed frames are also returned to the frames mimicking the behaviour of the get_frames().

Parameters:
  • frames (int or list or touple or str) – Specified frame (int), or frames (list), or range (touple), or all/everything (str). (default=’all’)
  • override (bool) – If True, an output already storred in analysis_output can be override. (default=False)
  • swap_atoms (dict, optional) – If this kwarg is passed with an appropriate dictionary a pywindow.molecular.MolecularSystem.swap_atom_keys() will be applied to the extracted frame.
  • forcefield (str, optional) – If this kwarg is passed with appropriate forcefield keyword a pywindow.molecular.MolecularSystem.decipher_atom_keys() will be applied to the extracted frame.
  • modular (bool, optional) – If this kwarg is passed a pywindow.molecular.MolecularSystem.make_modular() will be applied to the extracted frame. (default=False)
  • rebuild (bool, optional) – If this kwarg is passed a rebuild=True is passed to pywindow.molecular.MolecularSystem.make_modular() that will be applied to the extracted frame. (default=False)
  • ncpus (int, optional) – If ncpus > 1, then the analysis is performed in parallel for the specified number of parallel jobs. Otherwise, it runs in serial. (default=1)
Returns:

None – The function returns None, the analysis output is returned to analysis_output dictionary.

Return type:

NoneType

get_frames(frames='all', override=False, **kwargs)[source]

Extract frames from the trajectory file.

Depending on the passed parameters a frame, a list of particular frames, a range of frames (from, to), or all frames can be extracted with this function.

Parameters:
  • frames (int or list or touple or str) – Specified frame (int), or frames (list), or range (touple), or all/everything (str). (default=`all`)
  • override (bool) – If True, a frame already storred in frames can be override. (default=False)
  • extract_data (bool, optional) – If False, a frame is returned as a str block as in the trajectory file. Ohterwise, it is extracted and returned as pywindow.molecular.MolecularSystem. (default=True)
  • swap_atoms (dict, optional) – If this kwarg is passed with an appropriate dictionary a pywindow.molecular.MolecularSystem.swap_atom_keys() will be applied to the extracted frame.
  • forcefield (str, optional) – If this kwarg is passed with appropriate forcefield keyword a pywindow.molecular.MolecularSystem.decipher_atom_keys() will be applied to the extracted frame.
Returns:

save_analysis(filepath=None, **kwargs)[source]

Dump the content of analysis_output as JSON dictionary.

Parameters:filepath (str) – The filepath for the JSON file.
Returns:None
Return type:NoneType
save_frames(frames, filepath=None, filetype='pdb', **kwargs)[source]
class pywindow.trajectory.PDB(filepath)[source]

Bases: object

analysis(frames='all', ncpus=1, override=False, **kwargs)[source]

Perform structural analysis on a frame/ set of frames.

Depending on the passed parameters a frame, a list of particular frames, a range of frames (from, to), or all frames can be analysed with this function.

The analysis is performed on each frame and each discrete molecule in that frame separately. The steps are as follows:

  1. A frame is extracted and returned as a MolecularSystem.
  2. If swap_atoms is set the atom ids are swapped.
  3. If forcefield is set the atom ids are deciphered.
  4. If rebuild is set the molecules in the system are rebuild.
  5. Each discrete molecule is extracted as Molecule
  6. Each molecule is analysed with Molecule.full_analysis()
  7. Analysis output populates the analysis_output dictionary.

As the analysis of trajectories often have to be unique, many options are conditional.

A side effect of this function is that the analysed frames are also returned to the frames mimicking the behaviour of the get_frames().

Parameters:
  • frames (int or list or touple or str) – Specified frame (int), or frames (list), or range (touple), or all/everything (str). (default=’all’)
  • override (bool) – If True, an output already storred in analysis_output can be override. (default=False)
  • swap_atoms (dict, optional) – If this kwarg is passed with an appropriate dictionary a pywindow.molecular.MolecularSystem.swap_atom_keys() will be applied to the extracted frame.
  • forcefield (str, optional) – If this kwarg is passed with appropriate forcefield keyword a pywindow.molecular.MolecularSystem.decipher_atom_keys() will be applied to the extracted frame.
  • modular (bool, optional) – If this kwarg is passed a pywindow.molecular.MolecularSystem.make_modular() will be applied to the extracted frame. (default=False)
  • rebuild (bool, optional) – If this kwarg is passed a rebuild=True is passed to pywindow.molecular.MolecularSystem.make_modular() that will be applied to the extracted frame. (default=False)
  • ncpus (int, optional) – If ncpus > 1, then the analysis is performed in parallel for the specified number of parallel jobs. Otherwise, it runs in serial. (default=1)
Returns:

None – The function returns None, the analysis output is returned to analysis_output dictionary.

Return type:

NoneType

get_frames(frames='all', override=False, **kwargs)[source]

Extract frames from the trajectory file.

Depending on the passed parameters a frame, a list of particular frames, a range of frames (from, to), or all frames can be extracted with this function.

Parameters:
  • frames (int or list or touple or str) – Specified frame (int), or frames (list), or range (touple), or all/everything (str). (default=`all`)
  • override (bool) – If True, a frame already storred in frames can be override. (default=False)
  • extract_data (bool, optional) – If False, a frame is returned as a str block as in the trajectory file. Ohterwise, it is extracted and returned as pywindow.molecular.MolecularSystem. (default=True)
  • swap_atoms (dict, optional) – If this kwarg is passed with an appropriate dictionary a pywindow.molecular.MolecularSystem.swap_atom_keys() will be applied to the extracted frame.
  • forcefield (str, optional) – If this kwarg is passed with appropriate forcefield keyword a pywindow.molecular.MolecularSystem.decipher_atom_keys() will be applied to the extracted frame.
Returns:

save_analysis(filepath=None, **kwargs)[source]

Dump the content of analysis_output as JSON dictionary.

Parameters:filepath (str) – The filepath for the JSON file.
Returns:None
Return type:NoneType
save_frames(frames, filepath=None, filetype='pdb', **kwargs)[source]
class pywindow.trajectory.XYZ(filepath)[source]

Bases: object

A container for an XYZ type trajectory.

This function takes an XYZ trajectory file and maps it for the binary points in the file where each frame starts/ends. This way the process is fast, as it not require loading the trajectory into computer memory. When a frame is being extracted, it is only this frame that gets loaded to the memory.

Frames can be accessed individually and loaded as an unmodified string, returned as a pywindow.molecular.MolecularSystem (and analysed), dumped as PDB or XYZ or JSON (if dumped as a pywindow.molecular.MolecularSystem.system)

filepath

str – The filepath.

filename

str – The filename.

system_id

str – The system id inherited from the filename.

frames

dict – A dictionary that is populated, on the fly, with the extracted frames.

analysis_output

dict – A dictionary that is populated, on the fly, with the analysis output.

analysis(frames='all', ncpus=1, override=False, **kwargs)[source]

Perform structural analysis on a frame/ set of frames.

Depending on the passed parameters a frame, a list of particular frames, a range of frames (from, to), or all frames can be analysed with this function.

The analysis is performed on each frame and each discrete molecule in that frame separately. The steps are as follows:

  1. A frame is extracted and returned as a MolecularSystem.
  2. If swap_atoms is set the atom ids are swapped.
  3. If forcefield is set the atom ids are deciphered.
  4. If rebuild is set the molecules in the system are rebuild.
  5. Each discrete molecule is extracted as Molecule
  6. Each molecule is analysed with Molecule.full_analysis()
  7. Analysis output populates the analysis_output dictionary.

As the analysis of trajectories often have to be unique, many options are conditional.

A side effect of this function is that the analysed frames are also returned to the frames mimicking the behaviour of the get_frames().

Parameters:
  • frames (int or list or touple or str) – Specified frame (int), or frames (list), or range (touple), or all/everything (str). (default=’all’)
  • override (bool) – If True, an output already storred in analysis_output can be override. (default=False)
  • swap_atoms (dict, optional) – If this kwarg is passed with an appropriate dictionary a pywindow.molecular.MolecularSystem.swap_atom_keys() will be applied to the extracted frame.
  • forcefield (str, optional) – If this kwarg is passed with appropriate forcefield keyword a pywindow.molecular.MolecularSystem.decipher_atom_keys() will be applied to the extracted frame.
  • modular (bool, optional) – If this kwarg is passed a pywindow.molecular.MolecularSystem.make_modular() will be applied to the extracted frame. (default=False)
  • rebuild (bool, optional) – If this kwarg is passed a rebuild=True is passed to pywindow.molecular.MolecularSystem.make_modular() that will be applied to the extracted frame. (default=False)
  • ncpus (int, optional) – If ncpus > 1, then the analysis is performed in parallel for the specified number of parallel jobs. Otherwise, it runs in serial. (default=1)
Returns:

None – The function returns None, the analysis output is returned to analysis_output dictionary.

Return type:

NoneType

get_frames(frames='all', override=False, **kwargs)[source]

Extract frames from the trajectory file.

Depending on the passed parameters a frame, a list of particular frames, a range of frames (from, to), or all frames can be extracted with this function.

Parameters:
  • frames (int or list or touple or str) – Specified frame (int), or frames (list), or range (touple), or all/everything (str). (default=`all`)
  • override (bool) – If True, a frame already storred in frames can be override. (default=False)
  • extract_data (bool, optional) – If False, a frame is returned as a str block as in the trajectory file. Ohterwise, it is extracted and returned as pywindow.molecular.MolecularSystem. (default=True)
  • swap_atoms (dict, optional) – If this kwarg is passed with an appropriate dictionary a pywindow.molecular.MolecularSystem.swap_atom_keys() will be applied to the extracted frame.
  • forcefield (str, optional) – If this kwarg is passed with appropriate forcefield keyword a pywindow.molecular.MolecularSystem.decipher_atom_keys() will be applied to the extracted frame.
Returns:

save_analysis(filepath=None, **kwargs)[source]

Dump the content of analysis_output as JSON dictionary.

Parameters:filepath (str) – The filepath for the JSON file.
Returns:None
Return type:NoneType
save_frames(frames, filepath=None, filetype='pdb', **kwargs)[source]
pywindow.trajectory.make_supercell(system, matrix, supercell=[1, 1, 1])[source]

Return a supercell.

This functions takes the input unitcell and creates a supercell of it that is returned as a new pywindow.molecular.MolecularSystem.

Parameters:
  • system (pywindow.molecular.MolecularSystem.system) – The unit cell for creation of the supercell
  • matrix (numpy.array) – The unit cell parameters in form of a lattice.
  • supercell (list, optional) – A list that specifies the size of the supercell in the a, b and c direction. (default=[1, 1, 1])
Returns:

Returns the created supercell as a new MolecularSystem.

Return type:

pywindow.molecular.MolecularSystem

pywindow.utilities module

Module containing all general purpose functions shared by other modules.

This module is not intended for the direct use by a User. Therefore, I will only docstring functions if I see fit to do so.

LOG

11/07/18
Changed the way vector path is analysed. Now, the initial analysis is done with the geometrical formula for line-sphere intersection. Only the remaining vestors that do not intersect any van der Waals spheres are then analysed in the old way.
27/07/17
Fixed the cartesian coordinates -> fractional coordinates -> cartesian coordinates conversion related functions, creation of lattice array from unit cell parameters (triclinic system: so applicable to any) and conversion back to unit cell parameters. WORKS! inspiration from: http://www.ruppweb.org/Xray/tutorial/Coordinate%20system%20transformation.htm
26/07/17
Changed the way bonds are determined. Now, rather then fixed value a formula and covalent radii are used as explained in the Elemental_Radii spreadsheet (see tables module).

TO DO LIST

  • Fix and validate calculating shape descriptors: asphericity, acylindricity and the realtive shape anisotropy. (Not working at the moment)
  • In the find_windows() function, maybe change the way the EPS value for the DBSCAN() is estimates. Need to look how the distances change with the increase in size of the sampling sphere. (validate this with the MongoDB)
pywindow.utilities.acylidricity(S)[source]
pywindow.utilities.align_principal_ax(elements, coordinates)[source]
pywindow.utilities.angle_between_vectors(x, y)[source]

Calculate the angle between two vectors x and y.

pywindow.utilities.asphericity(S)[source]
pywindow.utilities.calc_acylidricity(elements, coordinates)[source]
pywindow.utilities.calc_asphericity(elements, coordinates)[source]
pywindow.utilities.calc_relative_shape_anisotropy(elements, coordinates)[source]
pywindow.utilities.calculate_pore_shape(elements, coordinates, adjust=1, increment=0.1, **kwargs)[source]

Return average diameter for a molecule.

pywindow.utilities.calculate_window_diameter(window, elements, coordinates, **kwargs)[source]
pywindow.utilities.cart2frac_all(coordinates, lattice_array)[source]

Convert all cartesian coordinates to fractional.

pywindow.utilities.cartisian_from_fractional(coordinate, lattice_array)[source]

Return cartesian coordinate from a fractional one.

pywindow.utilities.center_of_coor(coordinates)[source]

Return the centre of coordinates.

Parameters:coordinates (numpy.ndarray) – An array containing molecule’s coordinates.
Returns:An 1d array with coordinates of the centre of coordinates excluding elements’ masses.
Return type:numpy.ndarray
pywindow.utilities.center_of_mass(elements, coordinates)[source]

Return the centre of mass (COM).

Parameters:
  • elements (numpy.ndarray) – An array of all elements (type: str) in a molecule.
  • coordinates (numpy.ndarray) – An array containing molecule’s coordinates.
Returns:

An 1d array with coordinates of the centre of mass including elements’ masses.

Return type:

numpy.ndarray

pywindow.utilities.circumcircle(coordinates, atom_sets)[source]
pywindow.utilities.circumcircle_window(coordinates, atom_set)[source]
pywindow.utilities.compose_atom_list(*args)[source]

Return an atom list from elements and/or atom ids and coordinates.

An atom list is a special object that some pywindowfunctions uses. It is a nested list of lists with each individual list containing:

  1. [[element, coordinates (x, y, z)], …]
  2. [[element, atom key, coordinates (x, y, z)], …]

They work better for molecular re-building than two separate arrays for elements and coordinates do.

Parameters:
  • elements (numpy.ndarray) – An array of all elements (type: str) in a molecule.
  • coordinates (numpy.ndarray) – An array containing molecule’s coordinates.
  • atom_ids (numpy.ndarray, optional) – An array of all forcfield dependent atom keys (type:str) in a molecule.
Returns:

Version 1 or version 2 atom list depending on input parameters.

Return type:

list

Raises:

_FunctionError : :class:`Exception` – Raised when wrong number of parameters is passed to the function.

pywindow.utilities.correct_pore_diameter(com, *params)[source]

Return negative of a pore diameter. (optimisation function).

pywindow.utilities.create_supercell(system, supercell=[[-1, 1], [-1, 1], [-1, 1]])[source]

Create a supercell.

pywindow.utilities.decipher_atom_key(atom_key, forcefield)[source]

Return element for deciphered atom key.

This functions checks if the forcfield specified by user is supported and passes the atom key to the appropriate function for deciphering.

Parameters:
  • atom_key (str) – The atom key which is to be deciphered.
  • forcefield (str) – The forcefield to which the atom key belongs to.
Returns:

A string that is the periodic table element equvalent of forcefield atom key.

Return type:

str

pywindow.utilities.decompose_atom_list(atom_list)[source]

Return elements and/or atom ids and coordinates from an atom list.

Depending on input type of an atom list (version 1 or 2)

  1. [[element, coordinates (x, y, z)], …]
  2. [[element, atom key, coordinates (x, y, z)], …]

the function reverses what pywindow.utilities.compose_atom_list() do.

Parameters:atom_list (list) – A nested list of lists (version 1 or 2)
Returns:A touple of elements and coordinates arrays, or if input contained atom ideas, also atom ids array.
Return type:touple
pywindow.utilities.discrete_molecules(system, rebuild=None, tol=0.4)[source]

Decompose molecular system into individual discreet molecules.

Note

New formula for bonds: (26/07/17) The two atoms, x and y, are considered bonded if the distance between them, calculated with distance matrix, is within the ranges: .. :math:

Rcov(x) + Rcov(y) - t < R(x,y) < Rcov(x) + Rcov(y) + t

where Rcov is the covalent radius and the tolarenace (t) is set to 0.4 Angstrom.

pywindow.utilities.distance(a, b)[source]

Return the distance between two vectors (points) a and b.

Parameters:
  • a (numpy.ndarray) – First vector.
  • b (numpy.ndarray) – Second vector.
Returns:

A distance between two vectors (points).

Return type:

numpy.float64

pywindow.utilities.dlf_notation(atom_key)[source]

Return element for atom key using DL_F notation.

pywindow.utilities.find_average_diameter(elements, coordinates, adjust=1, increment=0.1, processes=None, **kwargs)[source]

Return average diameter for a molecule.

pywindow.utilities.find_windows(elements, coordinates, processes=None, mol_size=None, adjust=1, pore_opt=True, increment=1.0, **kwargs)[source]

Return windows diameters and center of masses for a molecule.

pywindow.utilities.find_windows_new(elements, coordinates, processes=None, mol_size=None, adjust=1, pore_opt=True, increment=1.0, **kwargs)[source]

Return windows diameters and center of masses for a molecule.

pywindow.utilities.frac2cart_all(frac_coordinates, lattice_array)[source]

Convert all fractional coordinates to cartesian.

pywindow.utilities.fractional_from_cartesian(coordinate, lattice_array)[source]

Return a fractional coordinate from a cartesian one.

pywindow.utilities.get_gyration_tensor(elements, coordinates)[source]

Return the gyration tensor of a molecule.

The gyration tensor should be invariant to the molecule’s position. The known formulas for the gyration tensor have the correction for the centre of mass of the molecule, therefore, the coordinates are first corrected for the centre of mass and essentially shifted to the origin.

Parameters:
  • elements (numpy.ndarray) – The array containing the molecule’s elemental data.
  • coordinates (numpy.ndarray) – The array containing the Cartesian coordinates of the molecule.
Returns:

The gyration tensor of a molecule invariant to the molecule’s position.

Return type:

numpy.ndarray

pywindow.utilities.get_inertia_tensor(elements, coordinates)[source]

Return the tensor of inertia a molecule.

Parameters:
  • elements (numpy.ndarray) – The array containing the molecule’s elemental data.
  • coordinates (numpy.ndarray) – The array containing the Cartesian coordinates of the molecule.
Returns:

The tensor of inertia of a molecule.

Return type:

numpy.ndarray

pywindow.utilities.get_tensor_eigenvalues(T, sort=False)[source]
pywindow.utilities.get_window_com(window, elements, coordinates, initial_com, **kwargs)[source]
pywindow.utilities.is_inside_polyhedron(point, polyhedron)[source]
pywindow.utilities.is_number(number)[source]

Return True if an object is a number - can be converted into a float.

Parameters:number (any) –
Returns:True if input is a float convertable (a number), False otherwise.
Return type:bool
pywindow.utilities.lattice_array_to_unit_cell(lattice_array)[source]

Return crystallographic param. from unit cell lattice matrix.

pywindow.utilities.max_dim(elements, coordinates)[source]

Return the maximum diameter of a molecule.

Parameters:
  • elements (numpy.ndarray) – An array of all elements (type: str) in a molecule.
  • coordinates (numpy.ndarray) – An array containing molecule’s coordinates.
pywindow.utilities.molecular_weight(elements)[source]

Return molecular weight of a molecule.

Parameters:elements (numpy.ndarray) – An array of all elements (type: str) in a molecule.
Returns:A molecular weight of a molecule.
Return type:numpy.float64
pywindow.utilities.normal_vector(origin, vectors)[source]

Return normal vector for two vectors with same origin.

pywindow.utilities.normalize_vector(vector)[source]

Normalize a vector.

A new vector is returned, the original vector is not modified.

Parameters:vector (np.array) – The vector to be normalized.
Returns:The normalized vector.
Return type:np.array
pywindow.utilities.opls_notation(atom_key)[source]

Return element for OPLS forcefield atom key.

pywindow.utilities.opt_pore_diameter(elements, coordinates, bounds=None, com=None, **kwargs)[source]

Return optimised pore diameter and it’s COM.

pywindow.utilities.optimise_xy(xy, *args)[source]

Return negative pore diameter for x and y coordinates optimisation.

pywindow.utilities.optimise_z(z, *args)[source]

Return pore diameter for coordinates optimisation in z direction.

pywindow.utilities.pore_diameter(elements, coordinates, com=None)[source]

Return pore diameter of a molecule.

pywindow.utilities.principal_axes(elements, coordinates)[source]
pywindow.utilities.relative_shape_anisotropy(S)[source]
pywindow.utilities.rotation_matrix_arbitrary_axis(angle, axis)[source]

Return a rotation matrix of angle radians about axis.

Parameters:
  • angle (int or float) – The size of the rotation in radians.
  • axis (numpy.array) – A 3 element aray which represents a vector. The vector is the axis about which the rotation is carried out.
Returns:

A 3x3 array representing a rotation matrix.

Return type:

numpy.array

pywindow.utilities.shift_com(elements, coordinates, com_adjust=array([ 0., 0., 0.]))[source]

Return coordinates translated by some vector.

Parameters:
  • elements (numpy.ndarray) – An array of all elements (type: str) in a molecule.
  • coordinates (numpy.ndarray) – An array containing molecule’s coordinates.
  • com_adjust (numpy.ndarray (default = [0, 0, 0])) –
Returns:

Translated array of molecule’s coordinates.

Return type:

numpy.ndarray

pywindow.utilities.sphere_volume(sphere_radius)[source]

Return volume of a sphere.

pywindow.utilities.to_list(obj)[source]
pywindow.utilities.unique(input_list)[source]

Return a list of unique items (similar to set functionality).

Parameters:input_list (list) – A list containg some items that can occur more than once.
Returns:A list with only unique occurances of an item.
Return type:list
pywindow.utilities.unit_cell_to_lattice_array(cryst)[source]

Return parallelpiped unit cell lattice matrix.

pywindow.utilities.vector_analysis(vector, coordinates, elements_vdw, increment=1.0)[source]

Analyse a sampling vector’s path for window analysis purpose.

pywindow.utilities.vector_analysis_pore_shape(vector, coordinates, elements_vdw)[source]
pywindow.utilities.vector_analysis_reversed(vector, coordinates, elements_vdw)[source]
pywindow.utilities.vector_preanalysis(vector, coordinates, elements_vdw, increment=1.0)[source]
pywindow.utilities.volume_from_cell_parameters(cryst)[source]

Return unit cell’s volume from crystallographic parameters.

pywindow.utilities.volume_from_lattice_array(lattice_array)[source]

Return unit cell’s volume from lattice matrix.

pywindow.utilities.window_analysis(window, elements, coordinates, elements_vdw, increment2=0.1, z_bounds=[None, None], lb_z=True, z_second_mini=False, **kwargs)[source]

Return window diameter and window’s centre.

Parameters:
  • widnow (list) –
  • elements (numpy.array) –
  • coordinates (numpy.array) –
  • elements_vdw (numpy.array) –
  • step (float) –
pywindow.utilities.window_shape(window, elements, coordinates, increment2=0.1, z_bounds=[None, None], lb_z=True, z_second_mini=False, **kwargs)[source]

Return window diameter and window’s centre.

Parameters:
  • widnow (list) –
  • elements (numpy.array) –
  • coordinates (numpy.array) –
  • elements_vdw (numpy.array) –
  • step (float) –

Module contents