ForceBalance API
1.3
Automated optimization of force fields and empirical potentials
|
Lee-Ping's general file format conversion class. More...
Public Member Functions | |
def | __init__ (self, fnm=None, ftype=None, top=None, ttype=None, kwargs) |
Create a Molecule object. More... | |
def | __len__ (self) |
Return the number of frames in the trajectory. More... | |
def | __getattr__ (self, key) |
Whenever we try to get a class attribute, it first tries to get the attribute from the Data dictionary. More... | |
def | __setattr__ (self, key, value) |
Whenever we try to get a class attribute, it first tries to get the attribute from the Data dictionary. More... | |
def | __deepcopy__ (self, memo) |
Custom deepcopy method because Python 3.6 appears to have changed its behavior. More... | |
def | __getitem__ (self, key) |
The Molecule class has list-like behavior, so we can get slices of it. More... | |
def | __delitem__ (self, key) |
Similarly, in order to delete a frame, we simply perform item deletion on framewise variables. More... | |
def | __iter__ (self) |
List-like behavior for looping over trajectories. More... | |
def | __add__ (self, other) |
Add method for Molecule objects. More... | |
def | __iadd__ (self, other) |
Add method for Molecule objects. More... | |
def | repair (self, key, klast) |
Attempt to repair trivial issues that would otherwise break the object. More... | |
def | reorder_according_to (self, other) |
Reorder atoms according to some other Molecule object. More... | |
def | reorder_indices (self, other) |
Return the indices that would reorder atoms according to some other Molecule object. More... | |
def | append (self, other) |
def | require (self, args) |
def | write (self, fnm=None, ftype=None, append=False, selection=None, kwargs) |
def | center_of_mass (self) |
def | radius_of_gyration (self) |
def | rigid_water (self) |
If one atom is oxygen and the next two are hydrogen, make the water molecule rigid. More... | |
def | load_frames (self, fnm, ftype=None, kwargs) |
def | edit_qcrems (self, in_dict, subcalc=None) |
Edit Q-Chem rem variables with a dictionary. More... | |
def | add_quantum (self, other) |
def | add_virtual_site (self, idx, kwargs) |
Add a virtual site to the system. More... | |
def | replace_peratom (self, key, orig, want) |
Replace all of the data for a certain attribute in the system from orig to want. More... | |
def | replace_peratom_conditional (self, key1, cond, key2, orig, want) |
Replace all of the data for a attribute key2 from orig to want, contingent on key1 being equal to cond. More... | |
def | atom_select (self, atomslice, build_topology=True) |
Return a copy of the object with certain atoms selected. More... | |
def | atom_stack (self, other) |
Return a copy of the object with another molecule object appended. More... | |
def | align_by_moments (self) |
Align molecules using the moment of inertia. More... | |
def | get_populations (self) |
Return a cloned molecule object but with X-coordinates set to Mulliken charges and Y-coordinates set to Mulliken spins. More... | |
def | load_popxyz (self, fnm) |
Given a charge-spin xyz file, load the charges (x-coordinate) and spins (y-coordinate) into internal arrays. More... | |
def | align (self, smooth=False, center=True, center_mass=False, atom_select=None) |
Align molecules. More... | |
def | center (self, center_mass=False) |
Move geometric center to the origin. More... | |
def | build_bonds (self) |
Build the bond connectivity graph. More... | |
def | build_topology (self, force_bonds=True, kwargs) |
Create self.topology and self.molecules; these are graph representations of the individual molecules (fragments) contained in the Molecule object. More... | |
def | distance_matrix (self, pbc=True) |
Obtain distance matrix between all pairs of atoms. More... | |
def | distance_displacement (self) |
Obtain distance matrix and displacement vectors between all pairs of atoms. More... | |
def | rotate_bond (self, frame, aj, ak, increment=15) |
Return a new Molecule object containing the selected frame plus a number of frames where the selected dihedral angle is rotated in steps of 'increment' given in degrees. More... | |
def | find_clashes (self, thre=0.0, pbc=True, groups=None) |
Obtain a list of atoms that 'clash' (i.e. More... | |
def | rotate_check_clash (self, frame, rotate_index, thresh_hyd=1.4, thresh_hvy=1.8, printLevel=1) |
Return a new Molecule object containing the selected frame plus a number of frames where the selected dihedral angle is rotated in steps of 'increment' given in degrees. More... | |
def | find_angles (self) |
Return a list of 3-tuples corresponding to all of the angles in the system. More... | |
def | find_dihedrals (self) |
Return a list of 4-tuples corresponding to all of the dihedral angles in the system. More... | |
def | measure_distances (self, i, j) |
def | measure_angles (self, i, j, k) |
def | measure_dihedrals (self, i, j, k, l) |
Return a series of dihedral angles, given four atom indices numbered from zero. More... | |
def | find_rings (self, max_size=6) |
Return a list of rings in the molecule. More... | |
def | order_by_connectivity (self, m, i, currList, max_min_path) |
Recursive function that orders atoms based on connectivity. More... | |
def | aliphatic_hydrogens (self) |
def | all_pairwise_rmsd (self) |
Find pairwise RMSD (super slow, not like the one in MSMBuilder.) More... | |
def | pathwise_rmsd (self, align=True) |
Find RMSD between frames along path. More... | |
def | ref_rmsd (self, i, align=True) |
Find RMSD to a reference frame. More... | |
def | align_center (self) |
def | openmm_positions (self) |
Returns the Cartesian coordinates in the Molecule object in a list of OpenMM-compatible positions, so it is possible to type simulation.context.setPositions(Mol.openmm_positions()[0]) or something like that. More... | |
def | openmm_boxes (self) |
Returns the periodic box vectors in the Molecule object in a list of OpenMM-compatible boxes, so it is possible to type simulation.context.setPeriodicBoxVectors(Mol.openmm_boxes()[0]) or something like that. More... | |
def | split (self, fnm=None, ftype=None, method="chunks", num=None) |
Split the molecule object into a number of separate files (chunks), either by specifying the number of frames per chunk or the number of chunks. More... | |
def | without (self, args) |
Return a copy of the Molecule object but with certain fields deleted if they exist. More... | |
def | read_comm_charge_mult (self, verbose=False) |
Set charge and multiplicity from reading the comment line, formatted in a specific way. More... | |
def | read_qcschema (self, schema, kwargs) |
def | read_xyz (self, fnm, kwargs) |
.xyz files can be TINKER formatted which is why we have the try/except here. More... | |
def | read_xyz0 (self, fnm, kwargs) |
Parse a .xyz file which contains several xyz coordinates, and return their elements. More... | |
def | read_mdcrd (self, fnm, kwargs) |
Parse an AMBER .mdcrd file. More... | |
def | read_inpcrd (self, fnm, kwargs) |
Parse an AMBER .inpcrd or .rst file. More... | |
def | read_qdata (self, fnm, kwargs) |
def | read_mol2 (self, fnm, kwargs) |
def | read_dcd (self, fnm, kwargs) |
def | read_com (self, fnm, kwargs) |
Parse a Gaussian .com file and return a SINGLE-ELEMENT list of xyz coordinates (no multiple file support) More... | |
def | read_arc (self, fnm, kwargs) |
Read a TINKER .arc file. More... | |
def | read_gro (self, fnm, kwargs) |
Read a GROMACS .gro file. More... | |
def | read_charmm (self, fnm, kwargs) |
Read a CHARMM .cor (or .crd) file. More... | |
def | read_qcin (self, fnm, kwargs) |
Read a Q-Chem input file. More... | |
def | read_pdb (self, fnm, kwargs) |
Loads a PDB and returns a dictionary containing its data. More... | |
def | read_qcesp (self, fnm, kwargs) |
def | read_qcout (self, fnm, errok=None, kwargs) |
Q-Chem output file reader, adapted for our parser. More... | |
def | write_qcin (self, selection, kwargs) |
def | write_xyz (self, selection, kwargs) |
def | get_reaxff_atom_types (self) |
Return a list of element names which maps the LAMMPS atom types to the ReaxFF elements. More... | |
def | write_lammps_data (self, selection, kwargs) |
Write the first frame of the selection to a LAMMPS data file for the purpose of automatically initializing a LAMMPS simulation. More... | |
def | write_molproq (self, selection, kwargs) |
def | write_mdcrd (self, selection, kwargs) |
def | write_inpcrd (self, selection, sn=None, kwargs) |
def | write_arc (self, selection, kwargs) |
def | write_gro (self, selection, kwargs) |
def | write_dcd (self, selection, kwargs) |
def | write_pdb (self, selection, kwargs) |
def | write_qdata (self, selection, kwargs) |
Text quantum data format. More... | |
def | require_resid (self) |
def | require_resname (self) |
def | require_boxes (self) |
Public Attributes | |
Read_Tab | |
The table of file readers. More... | |
Write_Tab | |
The table of file writers. More... | |
Funnel | |
A funnel dictionary that takes redundant file types and maps them down to a few. More... | |
positive_resid | |
Creates entries like 'gromacs' : 'gromacs' and 'xyz' : 'xyz' in the Funnel. More... | |
built_bonds | |
top_settings | |
Topology settings. More... | |
Data | |
comms | |
Read in stuff if we passed in a file name, otherwise return an empty instance. More... | |
fout | |
I needed to add in this line because the DCD writer requires the file name, but the other methods don't. More... | |
qm_mulliken_charges | |
qm_mulliken_spins | |
topology | |
molecules | |
charge | |
mult | |
xyzs | |
atomname | |
resid | |
resname | |
boxes | |
Lee-Ping's general file format conversion class.
The purpose of this class is to read and write chemical file formats in a way that is convenient for research. There are highly general file format converters out there (e.g. catdcd, openbabel) but I find that writing my own class can be very helpful for specific purposes. Here are some things this class can do:
Special variables: These variables cannot be set manually because there is a special method associated with getting them.
na = The number of atoms. You'll get this if you use MyMol.na or MyMol['na']. ns = The number of snapshots. You'll get this if you use MyMol.ns or MyMol['ns'].
Unit system: Angstroms.
Definition at line 1182 of file molecule.py.
def src.molecule.Molecule.__init__ | ( | self, | |
fnm = None , |
|||
ftype = None , |
|||
top = None , |
|||
ttype = None , |
|||
kwargs | |||
) |
Create a Molecule object.
fnm : str, optional File name to create the Molecule object from. If provided, the file will be parsed and used to fill in the fields such as elem (elements), xyzs (coordinates) and so on. If ftype is not provided, will automatically try to determine file type from file extension. If not provided, will create an empty object. ftype : str, optional File type, corresponding to an entry in the internal table of known file types. Provide this if you have a nonstandard file extension or if you wish to force to invoke a particular parser. top : str, optional "Topology" file name. If provided, will build the Molecule using this file name, then "fnm" will load frames instead. ttype : str, optional "Typology" file type. build_topology : bool, optional Build the molecular topology consisting of: topology (overall connectivity graph), molecules (list of connected subgraphs), bonds (if not explicitly read in), default True toppbc : bool, optional Use periodic boundary conditions when building the molecular topology, default False The build_topology code will attempt to determine this intelligently. topframe : int, optional Provide a frame number for building the molecular topology, default first frame Fac : float, optional Multiplicative factor to covalent radii criterion for deciding whether two atoms are bonded Default value of 1.2 is reasonable, 1.4 will produce lots of bonds positive_resid : bool, optional If provided, enforce all positive resIDs.
Definition at line 1219 of file molecule.py.
def src.molecule.Molecule.__add__ | ( | self, | |
other | |||
) |
Add method for Molecule objects.
Definition at line 1496 of file molecule.py.
def src.molecule.Molecule.__deepcopy__ | ( | self, | |
memo | |||
) |
Custom deepcopy method because Python 3.6 appears to have changed its behavior.
Definition at line 1398 of file molecule.py.
def src.molecule.Molecule.__delitem__ | ( | self, | |
key | |||
) |
Similarly, in order to delete a frame, we simply perform item deletion on framewise variables.
Definition at line 1476 of file molecule.py.
def src.molecule.Molecule.__getattr__ | ( | self, | |
key | |||
) |
Whenever we try to get a class attribute, it first tries to get the attribute from the Data dictionary.
Definition at line 1347 of file molecule.py.
def src.molecule.Molecule.__getitem__ | ( | self, | |
key | |||
) |
The Molecule class has list-like behavior, so we can get slices of it.
If we say MyMolecule[0:10], then we'll return a copy of MyMolecule with frames 0 through 9.
Definition at line 1453 of file molecule.py.
def src.molecule.Molecule.__iadd__ | ( | self, | |
other | |||
) |
Add method for Molecule objects.
Definition at line 1543 of file molecule.py.
def src.molecule.Molecule.__iter__ | ( | self | ) |
List-like behavior for looping over trajectories.
Note that these values are returned by reference. Note that this is intended to be more efficient than getitem, so when we loop over a trajectory, it's best to go "for m in M" instead of "for i in range(len(M)): m = M[i]"
Definition at line 1485 of file molecule.py.
def src.molecule.Molecule.__len__ | ( | self | ) |
Return the number of frames in the trajectory.
Definition at line 1331 of file molecule.py.
def src.molecule.Molecule.__setattr__ | ( | self, | |
key, | |||
value | |||
) |
Whenever we try to get a class attribute, it first tries to get the attribute from the Data dictionary.
Definition at line 1386 of file molecule.py.
def src.molecule.Molecule.add_quantum | ( | self, | |
other | |||
) |
Definition at line 1813 of file molecule.py.
def src.molecule.Molecule.add_virtual_site | ( | self, | |
idx, | |||
kwargs | |||
) |
Add a virtual site to the system.
This does NOT set the position of the virtual site; it sits at the origin.
Definition at line 1826 of file molecule.py.
def src.molecule.Molecule.align | ( | self, | |
smooth = False , |
|||
center = True , |
|||
center_mass = False , |
|||
atom_select = None |
|||
) |
Align molecules.
Has the option to create smooth trajectories (align each frame to the previous one) or to align each frame to the first one.
Also has the option to remove the center of mass.
Provide a list of atom indices to align along selected atoms.
Definition at line 1994 of file molecule.py.
def src.molecule.Molecule.align_by_moments | ( | self | ) |
Align molecules using the moment of inertia.
Departs from MSMBuilder convention of using arithmetic mean for mass.
Definition at line 1952 of file molecule.py.
def src.molecule.Molecule.align_center | ( | self | ) |
def src.molecule.Molecule.aliphatic_hydrogens | ( | self | ) |
Definition at line 2791 of file molecule.py.
def src.molecule.Molecule.all_pairwise_rmsd | ( | self | ) |
Find pairwise RMSD (super slow, not like the one in MSMBuilder.)
Definition at line 2803 of file molecule.py.
def src.molecule.Molecule.append | ( | self, | |
other | |||
) |
Definition at line 1652 of file molecule.py.
def src.molecule.Molecule.atom_select | ( | self, | |
atomslice, | |||
build_topology = True |
|||
) |
Return a copy of the object with certain atoms selected.
Takes an integer, list or array as argument.
Definition at line 1868 of file molecule.py.
def src.molecule.Molecule.atom_stack | ( | self, | |
other | |||
) |
Return a copy of the object with another molecule object appended.
WARNING: This function may invalidate stuff like QM energies.
Definition at line 1905 of file molecule.py.
def src.molecule.Molecule.build_bonds | ( | self | ) |
Build the bond connectivity graph.
Definition at line 2034 of file molecule.py.
def src.molecule.Molecule.build_topology | ( | self, | |
force_bonds = True , |
|||
kwargs | |||
) |
Create self.topology and self.molecules; these are graph representations of the individual molecules (fragments) contained in the Molecule object.
force_bonds : bool Build the bonds from interatomic distances. If the user calls build_topology from outside, assume this is the default behavior. If creating a Molecule object using init, do not force the building of bonds by default (only build bonds if not read from file.) topframe : int, optional Provide the frame number used for reading the bonds. If not provided, this will be taken from the top_settings field. If provided, this will take priority and write the value into top_settings.
Definition at line 2208 of file molecule.py.
def src.molecule.Molecule.center | ( | self, | |
center_mass = False |
|||
) |
Move geometric center to the origin.
Definition at line 2023 of file molecule.py.
def src.molecule.Molecule.center_of_mass | ( | self | ) |
Definition at line 1723 of file molecule.py.
def src.molecule.Molecule.distance_displacement | ( | self | ) |
Obtain distance matrix and displacement vectors between all pairs of atoms.
Definition at line 2254 of file molecule.py.
def src.molecule.Molecule.distance_matrix | ( | self, | |
pbc = True |
|||
) |
Obtain distance matrix between all pairs of atoms.
Definition at line 2242 of file molecule.py.
def src.molecule.Molecule.edit_qcrems | ( | self, | |
in_dict, | |||
subcalc = None |
|||
) |
Edit Q-Chem rem variables with a dictionary.
Pass a value of None to delete a rem variable.
Definition at line 1798 of file molecule.py.
def src.molecule.Molecule.find_angles | ( | self | ) |
Return a list of 3-tuples corresponding to all of the angles in the system.
Verified for lysine and tryptophan dipeptide when comparing to TINKER's analyze program.
Definition at line 2502 of file molecule.py.
def src.molecule.Molecule.find_clashes | ( | self, | |
thre = 0.0 , |
|||
pbc = True , |
|||
groups = None |
|||
) |
Obtain a list of atoms that 'clash' (i.e.
are more than 3 bonds apart and are closer than the provided threshold.)
thre : float Create a sorted-list of all non-bonded atom pairs with distance below this threshold pbc : bool Whether to use PBC when computing interatomic distances filt : tuple (group1, group2) If not None, only compute distances between atoms in 'group1' and atoms in 'group2'. For example this could be [gAtoms, oAtoms] from rotate_bond
minPair_frames : list of tuples The closest pair of non-bonded atoms in each frame minDist_frames : list of tuples The distance between the closest pair of non-bonded atoms in each frame clashPairs_frames : list of lists Pairs of non-bonded atoms with distance below the threshold in each frame clashDists_frames : list of lists Distances between pairs of non-bonded atoms below the threshold in each frame
Definition at line 2379 of file molecule.py.
def src.molecule.Molecule.find_dihedrals | ( | self | ) |
Return a list of 4-tuples corresponding to all of the dihedral angles in the system.
Verified for alanine and tryptophan dipeptide when comparing to TINKER's analyze program.
Definition at line 2529 of file molecule.py.
def src.molecule.Molecule.find_rings | ( | self, | |
max_size = 6 |
|||
) |
Return a list of rings in the molecule.
Tested on a DNA base pair and C60. Warning: Using large max_size for rings (e.g. for finding the macrocycle in porphyrin) could lead to some undefined behavior.
max_size : int The maximum ring size. If a large ring contains smaller sub-rings, they are all mapped into one.
list A list of rings sorted in ascending order of the smallest atom number. The atoms in each ring are ordered starting from the lowest number and going along the ring, with the second atom being the lower of the two possible choices.
Definition at line 2621 of file molecule.py.
def src.molecule.Molecule.get_populations | ( | self | ) |
Return a cloned molecule object but with X-coordinates set to Mulliken charges and Y-coordinates set to Mulliken spins.
Definition at line 1966 of file molecule.py.
def src.molecule.Molecule.get_reaxff_atom_types | ( | self | ) |
Return a list of element names which maps the LAMMPS atom types to the ReaxFF elements.
Definition at line 4229 of file molecule.py.
def src.molecule.Molecule.load_frames | ( | self, | |
fnm, | |||
ftype = None , |
|||
kwargs | |||
) |
Definition at line 1769 of file molecule.py.
def src.molecule.Molecule.load_popxyz | ( | self, | |
fnm | |||
) |
Given a charge-spin xyz file, load the charges (x-coordinate) and spins (y-coordinate) into internal arrays.
Definition at line 1977 of file molecule.py.
def src.molecule.Molecule.measure_angles | ( | self, | |
i, | |||
j, | |||
k | |||
) |
Definition at line 2560 of file molecule.py.
def src.molecule.Molecule.measure_dihedrals | ( | self, | |
i, | |||
j, | |||
k, | |||
l | |||
) |
Return a series of dihedral angles, given four atom indices numbered from zero.
Definition at line 2575 of file molecule.py.
def src.molecule.Molecule.measure_distances | ( | self, | |
i, | |||
j | |||
) |
Definition at line 2551 of file molecule.py.
def src.molecule.Molecule.openmm_boxes | ( | self | ) |
Returns the periodic box vectors in the Molecule object in a list of OpenMM-compatible boxes, so it is possible to type simulation.context.setPeriodicBoxVectors(Mol.openmm_boxes()[0]) or something like that.
Definition at line 2880 of file molecule.py.
def src.molecule.Molecule.openmm_positions | ( | self | ) |
Returns the Cartesian coordinates in the Molecule object in a list of OpenMM-compatible positions, so it is possible to type simulation.context.setPositions(Mol.openmm_positions()[0]) or something like that.
Definition at line 2863 of file molecule.py.
def src.molecule.Molecule.order_by_connectivity | ( | self, | |
m, | |||
i, | |||
currList, | |||
max_min_path | |||
) |
Recursive function that orders atoms based on connectivity.
To call the function, pass in the topology object (e.g. M.molecules[0]) and the index of the atom assigned to be "first".
The neighbors of "i" are placed in decreasing order of the maximum length of the shortest path to all other atoms - i.e. an atom closer to the "edge" of the molecule has a longer shortest-path to atoms on the other "edge", and they should be added first.
Snippet:
from forcebalance.molecule import * import numpy as np M = Molecule('movProt.xyz', Fac=1.25)
firstAtom = [m.L()[np.argmax([PeriodicTable[M.elem[i]] for i in m.L()])] for m in M.molecules]
firstAtom[0] = 40
newOrder = [] for i in range(len(M.molecules)): newOrder += M.order_by_connectivity(M.molecules[i], firstAtom[i], [], None)
newM = M.atom_select(newOrder) newM.write('reordered1.xyz')
m : topology object (contained in Molecule) i : integer Index of the atom to be added first (if called recursively, the index of subsequent atoms to be added) currList : list Current list of new atom orderings max_min_path : dict Dictionary mapping
Definition at line 2758 of file molecule.py.
def src.molecule.Molecule.pathwise_rmsd | ( | self, | |
align = True |
|||
) |
Find RMSD between frames along path.
Definition at line 2821 of file molecule.py.
def src.molecule.Molecule.radius_of_gyration | ( | self | ) |
def src.molecule.Molecule.read_arc | ( | self, | |
fnm, | |||
kwargs | |||
) |
Read a TINKER .arc file.
[in] | fnm | The input file name |
Definition at line 3300 of file molecule.py.
def src.molecule.Molecule.read_charmm | ( | self, | |
fnm, | |||
kwargs | |||
) |
Read a CHARMM .cor (or .crd) file.
Definition at line 3455 of file molecule.py.
def src.molecule.Molecule.read_com | ( | self, | |
fnm, | |||
kwargs | |||
) |
Parse a Gaussian .com file and return a SINGLE-ELEMENT list of xyz coordinates (no multiple file support)
[in] | fnm | The input file name |
Definition at line 3255 of file molecule.py.
def src.molecule.Molecule.read_comm_charge_mult | ( | self, | |
verbose = False |
|||
) |
Set charge and multiplicity from reading the comment line, formatted in a specific way.
Definition at line 2925 of file molecule.py.
def src.molecule.Molecule.read_dcd | ( | self, | |
fnm, | |||
kwargs | |||
) |
def src.molecule.Molecule.read_gro | ( | self, | |
fnm, | |||
kwargs | |||
) |
Read a GROMACS .gro file.
Definition at line 3373 of file molecule.py.
def src.molecule.Molecule.read_inpcrd | ( | self, | |
fnm, | |||
kwargs | |||
) |
Parse an AMBER .inpcrd or .rst file.
[in] | fnm | The input file name |
Definition at line 3064 of file molecule.py.
def src.molecule.Molecule.read_mdcrd | ( | self, | |
fnm, | |||
kwargs | |||
) |
Parse an AMBER .mdcrd file.
This requires at least the number of atoms. This will FAIL for monatomic trajectories (but who the heck makes those?)
[in] | fnm | The input file name |
Definition at line 3031 of file molecule.py.
def src.molecule.Molecule.read_mol2 | ( | self, | |
fnm, | |||
kwargs | |||
) |
Definition at line 3161 of file molecule.py.
def src.molecule.Molecule.read_pdb | ( | self, | |
fnm, | |||
kwargs | |||
) |
Loads a PDB and returns a dictionary containing its data.
Definition at line 3649 of file molecule.py.
def src.molecule.Molecule.read_qcesp | ( | self, | |
fnm, | |||
kwargs | |||
) |
def src.molecule.Molecule.read_qcin | ( | self, | |
fnm, | |||
kwargs | |||
) |
Read a Q-Chem input file.
These files can be very complicated, and I can't write a completely general parser for them. It is important to keep our goal in mind:
1) The main goal is to convert a trajectory to Q-Chem input files with identical calculation settings.
2) When we print the Q-Chem file, we should preserve the line ordering of the 'rem' section, but also be able to add 'rem' options at the end.
3) We should accommodate the use case that the Q-Chem file may have follow-up calculations delimited by '@@'.
4) We can read in all of the xyz's as a trajectory, but only the Q-Chem settings belonging to the first xyz will be saved.
Definition at line 3524 of file molecule.py.
def src.molecule.Molecule.read_qcout | ( | self, | |
fnm, | |||
errok = None , |
|||
kwargs | |||
) |
Q-Chem output file reader, adapted for our parser.
Q-Chem output files are very flexible and there's no way I can account for all of them. Here's what I am able to account for:
A list of:
Calling with errok will proceed with reading file even if the specified error messages are encountered.
Note that each step in a geometry optimization counts as a frame.
As with all Q-Chem output files, note that successive calculations can have different numbers of atoms.
Definition at line 3778 of file molecule.py.
def src.molecule.Molecule.read_qcschema | ( | self, | |
schema, | |||
kwargs | |||
) |
Definition at line 2933 of file molecule.py.
def src.molecule.Molecule.read_qdata | ( | self, | |
fnm, | |||
kwargs | |||
) |
Definition at line 3125 of file molecule.py.
def src.molecule.Molecule.read_xyz | ( | self, | |
fnm, | |||
kwargs | |||
) |
.xyz files can be TINKER formatted which is why we have the try/except here.
Definition at line 2952 of file molecule.py.
def src.molecule.Molecule.read_xyz0 | ( | self, | |
fnm, | |||
kwargs | |||
) |
Parse a .xyz file which contains several xyz coordinates, and return their elements.
[in] | fnm | The input file name |
Definition at line 2967 of file molecule.py.
def src.molecule.Molecule.ref_rmsd | ( | self, | |
i, | |||
align = True |
|||
) |
Find RMSD to a reference frame.
Definition at line 2839 of file molecule.py.
def src.molecule.Molecule.reorder_according_to | ( | self, | |
other | |||
) |
Reorder atoms according to some other Molecule object.
This happens when we run a program like pdb2gmx or pdbxyz and it scrambles our atom ordering, forcing us to reorder the atoms for all frames in the current Molecule object.
Directions: (1) Load up the scrambled file as a new Molecule object. (2) Call this function: Original_Molecule.reorder_according_to(scrambled) (3) Save Original_Molecule to a new file name.
Definition at line 1621 of file molecule.py.
def src.molecule.Molecule.reorder_indices | ( | self, | |
other | |||
) |
Return the indices that would reorder atoms according to some other Molecule object.
This happens when we run a program like pdb2gmx or pdbxyz and it scrambles our atom ordering.
Directions: (1) Load up the scrambled file as a new Molecule object. (2) Call this function: Original_Molecule.reorder_indices(scrambled)
Definition at line 1648 of file molecule.py.
def src.molecule.Molecule.repair | ( | self, | |
key, | |||
klast | |||
) |
Attempt to repair trivial issues that would otherwise break the object.
Definition at line 1586 of file molecule.py.
def src.molecule.Molecule.replace_peratom | ( | self, | |
key, | |||
orig, | |||
want | |||
) |
Replace all of the data for a certain attribute in the system from orig to want.
Definition at line 1845 of file molecule.py.
def src.molecule.Molecule.replace_peratom_conditional | ( | self, | |
key1, | |||
cond, | |||
key2, | |||
orig, | |||
want | |||
) |
Replace all of the data for a attribute key2 from orig to want, contingent on key1 being equal to cond.
For instance: replace H1 with H2 if resname is SOL.
Definition at line 1857 of file molecule.py.
def src.molecule.Molecule.require | ( | self, | |
args | |||
) |
Definition at line 1656 of file molecule.py.
def src.molecule.Molecule.require_boxes | ( | self | ) |
def src.molecule.Molecule.require_resid | ( | self | ) |
Definition at line 4595 of file molecule.py.
def src.molecule.Molecule.require_resname | ( | self | ) |
Definition at line 4603 of file molecule.py.
def src.molecule.Molecule.rigid_water | ( | self | ) |
If one atom is oxygen and the next two are hydrogen, make the water molecule rigid.
Definition at line 1739 of file molecule.py.
def src.molecule.Molecule.rotate_bond | ( | self, | |
frame, | |||
aj, | |||
ak, | |||
increment = 15 |
|||
) |
Return a new Molecule object containing the selected frame plus a number of frames where the selected dihedral angle is rotated in steps of 'increment' given in degrees.
This function is designed to be called by Molecule.rotate_check_clash().
frame : int Structure number of the current Molecule to be rotated aj, ak : int Atom numbers of the bond to be rotated increment : float Degrees of the rotation increment
Molecule New Molecule object containing the rotated structures
Definition at line 2286 of file molecule.py.
def src.molecule.Molecule.rotate_check_clash | ( | self, | |
frame, | |||
rotate_index, | |||
thresh_hyd = 1.4 , |
|||
thresh_hvy = 1.8 , |
|||
printLevel = 1 |
|||
) |
Return a new Molecule object containing the selected frame plus a number of frames where the selected dihedral angle is rotated in steps of 'increment' given in degrees.
Additionally, check for if pairs of non-bonded atoms "clash" i.e. approach below the specified thresholds.
frame : int Structure number of the current Molecule to be rotated (if only a single frame, this number should be 0) rotate_index : tuple 4 atom indices (ai, aj, ak, al) representing the dihedral angle to be rotated. The first and last indices are used to measure the dihedral angle. thresh_hyd : float Clash threshold for hydrogen atoms. Reasonable values are in between 1.3 and 2.0. thresh_hvy : float Clash threshold for heavy atoms. Reasonable values are in between 1.7 and 2.5. printLevel: int Sets the amount of printout (larger = more printout)
Molecule New Molecule object containing the rotated structures Success : bool True if no clashes
Definition at line 2445 of file molecule.py.
def src.molecule.Molecule.split | ( | self, | |
fnm = None , |
|||
ftype = None , |
|||
method = "chunks" , |
|||
num = None |
|||
) |
Split the molecule object into a number of separate files (chunks), either by specifying the number of frames per chunk or the number of chunks.
Only relevant for "trajectories". The type of file may be specified; if they aren't specified then the original file type is used.
The output file names are [name].[numbers].[extension] where [name] can be specified by passing 'fnm' or taken from the object's 'fnm' attribute by default. [numbers] are integers ranging from the lowest to the highest chunk number, prepended by zeros.
If the number of chunks / frames is not specified, then one file is written for each frame.
Definition at line 2903 of file molecule.py.
def src.molecule.Molecule.without | ( | self, | |
args | |||
) |
Return a copy of the Molecule object but with certain fields deleted if they exist.
This is useful for when we'd like to add Molecule objects that don't share some fields that we know about.
Definition at line 2915 of file molecule.py.
def src.molecule.Molecule.write | ( | self, | |
fnm = None , |
|||
ftype = None , |
|||
append = False , |
|||
selection = None , |
|||
kwargs | |||
) |
def src.molecule.Molecule.write_arc | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_dcd | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_gro | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_inpcrd | ( | self, | |
selection, | |||
sn = None , |
|||
kwargs | |||
) |
def src.molecule.Molecule.write_lammps_data | ( | self, | |
selection, | |||
kwargs | |||
) |
Write the first frame of the selection to a LAMMPS data file for the purpose of automatically initializing a LAMMPS simulation.
This function makes several assumptions until further notice:
(1) We are interested in a ReaxFF simulation (2) Atom types will be generated from elements
Definition at line 4245 of file molecule.py.
def src.molecule.Molecule.write_mdcrd | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_molproq | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_pdb | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_qcin | ( | self, | |
selection, | |||
kwargs | |||
) |
def src.molecule.Molecule.write_qdata | ( | self, | |
selection, | |||
kwargs | |||
) |
Text quantum data format.
Definition at line 4574 of file molecule.py.
def src.molecule.Molecule.write_xyz | ( | self, | |
selection, | |||
kwargs | |||
) |
src.molecule.Molecule.atomname |
Definition at line 4458 of file molecule.py.
src.molecule.Molecule.boxes |
Definition at line 4658 of file molecule.py.
src.molecule.Molecule.built_bonds |
Definition at line 1281 of file molecule.py.
src.molecule.Molecule.charge |
Definition at line 2927 of file molecule.py.
src.molecule.Molecule.comms |
Read in stuff if we passed in a file name, otherwise return an empty instance.
Fill in comments.
Try to determine from the file name using the extension. Actually read the file. Set member variables. Create a list of comment lines if we don't already have them from reading the file.
Definition at line 1311 of file molecule.py.
src.molecule.Molecule.Data |
Definition at line 1293 of file molecule.py.
src.molecule.Molecule.fout |
I needed to add in this line because the DCD writer requires the file name, but the other methods don't.
Definition at line 1692 of file molecule.py.
src.molecule.Molecule.Funnel |
A funnel dictionary that takes redundant file types and maps them down to a few.
Definition at line 1264 of file molecule.py.
src.molecule.Molecule.molecules |
Definition at line 2235 of file molecule.py.
src.molecule.Molecule.mult |
Definition at line 2928 of file molecule.py.
src.molecule.Molecule.positive_resid |
Creates entries like 'gromacs' : 'gromacs' and 'xyz' : 'xyz' in the Funnel.
Definition at line 1280 of file molecule.py.
src.molecule.Molecule.qm_mulliken_charges |
Definition at line 1979 of file molecule.py.
src.molecule.Molecule.qm_mulliken_spins |
Definition at line 1980 of file molecule.py.
src.molecule.Molecule.Read_Tab |
The table of file readers.
Definition at line 1235 of file molecule.py.
src.molecule.Molecule.resid |
Definition at line 4599 of file molecule.py.
src.molecule.Molecule.resname |
Definition at line 4606 of file molecule.py.
src.molecule.Molecule.top_settings |
Topology settings.
Definition at line 1283 of file molecule.py.
src.molecule.Molecule.topology |
Definition at line 2233 of file molecule.py.
src.molecule.Molecule.Write_Tab |
The table of file writers.
Definition at line 1251 of file molecule.py.
src.molecule.Molecule.xyzs |
Definition at line 4338 of file molecule.py.