MM_Wrapper Submodule¶
janus.mm_wrapper.MMWrapper (class_type[, …]) |
A super class for all molecular mechanics wrappers |
janus.mm_wrapper.OpenMMWrapper ([sys_info, …]) |
A MM wrapper class that calls OpenMM to obtain molecular mechanics information. |
MMWrapper¶
-
class
janus.mm_wrapper.
MMWrapper
(class_type, sys_info=None, sys_info_format=None)¶ Bases:
abc.ABC
A super class for all molecular mechanics wrappers
Note
Since MMWrapper is a super class and has abstract methods the user cannot actually instantiate a MMWrapper object, but only its child objectsbuild_qm_param
()Function not implemented for MM wrappers compute_info
()Function implemented in individual child classes convert_input
()Function implemented in individual child classes convert_trajectory
()Function implemented in individual child classes get_energy_and_gradient
(traj[, geometry, …])Gets the energy and gradient from a MM computation get_geom_from_trajectory
()Function not implemented for MM wrappers get_main_charges
()Function implemented in individual child classes get_main_info
()Function implemented in individual child classes initialize
()Function implemented in individual child classes optimize_geometry
()Function not implemented for MM wrappers set_external_charges
()Function implemented in individual child classes set_qm_geometry
()Function not implemented for MM wrappers set_up_reporters
()Function implemented in individual child classes take_step
(force)Function implemented in individual child classes Attributes
au_bohr_to_kjmol_nm
kjmol_nm_to_au_bohr
kjmol_to_au
nm_to_angstrom
nm_to_bohr
-
au_bohr_to_kjmol_nm
= 49614.50168278496¶
-
build_qm_param
()¶ Function not implemented for MM wrappers
-
compute_info
()¶ Function implemented in individual child classes
-
convert_input
()¶ Function implemented in individual child classes
-
convert_trajectory
()¶ Function implemented in individual child classes
-
get_energy_and_gradient
(traj, geometry=None, include_coulomb='all', link_atoms=None, minimize=False, charges=None)¶ Gets the energy and gradient from a MM computation
Parameters: - traj (MDtraj trajectory object) –
- geometry (str) – A string containing geometry information as XYZ coordinates. Not applicable for MM programs
- include_coulomb (str) – whether to include coulombic interactions. ‘all’ (default) includes coulombic forces for all particles, ‘no_link’ excludes coulombic forces for link atoms, ‘only’ excludes all other forces for all atoms, ‘none’ excludes coulombic forces for all particles.
- link_atoms (list) – indices of link_atoms
- minimize (bool) – whether to return the geometry optimized energy
- charges (list) – charges and corresponding positions in angstroms as xyz coordinates
Returns: A dictionary with energy(‘energy’) and gradient(‘gradients’) information
Return type: dict
-
get_geom_from_trajectory
()¶ Function not implemented for MM wrappers
-
get_main_charges
()¶ Function implemented in individual child classes
-
get_main_info
()¶ Function implemented in individual child classes
-
initialize
()¶ Function implemented in individual child classes
-
kjmol_nm_to_au_bohr
= 2.0155397435886694e-05¶
-
kjmol_to_au
= 0.0003808798033989866¶
-
nm_to_angstrom
= 10.0¶
-
nm_to_bohr
= 0.052917999999¶
-
optimize_geometry
()¶ Function not implemented for MM wrappers
-
set_external_charges
()¶ Function implemented in individual child classes
-
set_qm_geometry
()¶ Function not implemented for MM wrappers
-
set_up_reporters
()¶ Function implemented in individual child classes
-
take_step
(force)¶ Function implemented in individual child classes
-
OpenMMWrapper¶
-
class
janus.mm_wrapper.
OpenMMWrapper
(sys_info=None, sys_info_format='pdb', mm_forcefield='amber99sb.xml', mm_water_forcefield='tip3p.xml', NVE_integrator='Verlet', NVT_integrator='Langevin', temp=300, step_size=1, fric_coeff=1, nonbondedCutoff=0.8, **kwargs)¶ Bases:
janus.mm_wrapper.mm_wrapper.MMWrapper
A MM wrapper class that calls OpenMM to obtain molecular mechanics information. Also can be used to take steps in a molecular dynamics simulation. Class inherits from MMWrapper.
Parameters: - sys_info (str) – A string with the filename or a list with multiple filenames that contain position and topology information. Default is None.
- sys_info_format (str) – A str describing what kind of input is contained in sys_info. Default is pdb. Possible values also include Amber and Gromacs.
- mm_forcefield (str) – The name of the forcefield to use, default is amber99sb.xml.
- mm_water_forcefield (str) – The name of the forcefield to use for water, default is tip3p.xml.
- NVE_integrator (str) – What type of integrator to use for a NVE ensemble, default is Verlet.
- NVT_integrator (str) – What type of integrator to use for a NVT ensemble, default is Langevin.
- temp (Int) – The temperature at which to run a simulation in kelvin, default is 300.
- step_size (float) – The step size to integrate system in femtoseconds, default is 1.
- fric_coeff (float) – friction coefficient to couple the system to heat bath in a NVT ensemble in inverse picoseconds, default is 1.
- nonbondedCutoff (float) – The cutoff distance for nonbonded interactions in nanometers, default is 1.
- **kwargs –
Other parameters for OpenMM, which include: - nonbondedMethod : method for nonbonded interactions, default is OM_app.NoCutoff - constraints : which bonds and angles implemented with constraints,
default is None- rigid_water : whether water is treated as rigid, default is True
- removeCMMotion : whether to include a CMMotionRemover, default is True
- ignoreExternalBonds : whether to ignore external bonds when matching residues to templates,
- default is True
- flexibleConstraints : whether to add parameters for constrained parameters,
- default is False
- hydrogenMass : the mass to use for hydrogen atoms bonded to heavy atoms,
- default is False
- residueTemplates : allows user to specify a template for a residue,
- default is empty dict {}
- switchDistance : the distance to turn on potential energy switching function for
- Lennard-Jones interactions. Default is None
- keywords for MD simulation parameters.
- For possible keywords consult the Molecular Dynamics section of the manual.
Each pair of key:value in the dictionary is given as a string. For more information about these pararmeters and other possible parameter values consult docs.openmm.org
-
compute_info
(topology, positions, include_coulomb='all', initialize=False, return_system=False, return_simulation=False, link_atoms=None, minimize=False)¶ Gets information about a system.
Parameters: - topology (OpenMM topology object) –
- positions (OpenMM Vec3 vector) – contains the positions of the system in nm
- include_coulomb (str) – whether to include coulombic interactions. ‘all’ (default) includes coulombic forces for all particles, ‘no_link’ excludes coulombic forces for link atoms, ‘only’ excludes all other forces for all atoms, ‘none’ excludes coulombic forces for all particles.
- initialize (bool) – Whether the main system is being initialized.
- return_system (bool) – True(default) to return OpenMM system object
- return_simulation (bool) – True(default) to return OpenMM simulation object
- link_atoms (list) – if included as a list with include_coulomb=’no_link’, specifies which atoms to remove coulombic forces from. Default is None.
- minimize (bool) – whether to minimize the energy of the system
Returns: - dict – A dictionary with state information
- OpenMM system object – OpenMM system object returned unless return_system=False
- OpenMM simulation object – OpenMM simulation object returned unless return_simulation=False
Examples
>>> system, simulation, state = compute_info(top, pos)
>>> state = compute_info(top, pos, return_simulation=False, return_system=False)
-
convert_input
()¶ Converts inputs to OpenMM readable topologies and positions. Currently supports pdb inputs as well as Amber and Gromacs input files.
-
convert_trajectory
(traj)¶ Converts an OpenMM trajectory to get topology and positions that are compatible with MDtraj
Parameters: traj (OpenMM trajectory object) – Returns: - list – positions in nm
- OpenMM topology object
Examples
>>> positions, topology = convert_trajectory(OpenMM_traj)
-
create_modeller
(atoms, keep_atoms=False)¶ Makes a OpenMM modeller object based on given geometry
Parameters: - atoms (list) – The subset of atom indices to either keep or delete
- keep_atoms (bool) – Whether to keep the atoms specified in the modeller or delete them. Default is false.
Returns: Return type: A OpenMM modeller object
Examples
>>> modeller = self.make_modeller() modeller = self.make_modeller(keep_qm=True)
-
create_new_residue_template
(topology)¶ Create a new OpeMM residue template when there is no matching residue and registers it into self.forcefield forcefield object.
Note
currently, if there is unmatched name, currently only checks original unmodified residue, N-terminus form, and C-terminus form. This may not be robust.
Parameters: topology (OpenMM topology object) – Examples
>>> create_new_residue_template(topology)
-
create_openmm_simulation
(openmm_system, topology, positions, integrator, return_integrator=False, seed=0)¶ Creates an OpenMM simulation object given an OpenMM system, topology, and positions
Parameters: - openmm_system (OpenMM system object) –
- topology (OpenMM topology object) –
- positions (OpenMM Vec3 vector) – contains the positions of the system in nm
- integrator (str) – What type of integrator to use. Currently support Langevin and Verlet.
- return_integrator (bool) – Whether to return the OpenMM integrator object
- seed (int) – Set a random seed number for the Langevin integrator. Default is 0, which means seed is randomized every time.
Returns: - OpenMM simulation object
- OpenMM integrator object – Returned unless return_integrator is False
Examples
>>> create_open_simulation(openmm_sys, top, pos) create_open_simulation(openmm_sys, pdb.topology. pdb.positions)
-
create_openmm_system
(topology, include_coulomb='all', link_atoms=None, initialize=False)¶ Calls OpenMM to create an OpenMM System object give a topology, forcefield, and other parameters specified in the instantiation parameters.
Note
Currently there is no support for customized forcefields
Parameters: - topology (OpenMM topology object) –
- include_coulomb (str) – whether to include coulombic interactions. ‘all’ (default) includes coulombic forces for all particles, ‘no_link’ excludes coulombic forces for link atoms, ‘only’ excludes all other forces for all atoms, ‘none’ excludes coulombic forces for all particles.
- link_atoms (list) – if included as a list with include_coulomb=’no_link’, specifies which atoms to remove coulombic forces from. Default is None.
- initialize (bool) – Whether the main system is being initialized.
Returns: Return type: OpenMM system object
Examples
>>> openmm_sys = create_openmm_system(topology) openmm_sys = create_openmm_system(top, include_coulomb='no_link', link_atoms=[0,1,2]) openmm_sys = create_openmm_system(top, include_coulomb='only') openmm_sys = create_openmm_system(top, initialize=True)
-
delete_atoms
(atoms)¶ Delete specified atoms from an OpenMM Modeller object
Parameters: - model (OpenMM Modeller object) –
- atoms (list) – which atom IDs (int) or atom names (str) to delete an OpenMM Modeller object
Examples
>>> delete_atoms(model, [0, 3, 5]) delete_atoms(model, ['Cl'])
-
get_main_charges
()¶ Gets the MM point charges for the system of interest
Returns: charges of system Return type: list Examples
>>> charges = get_main_charges()
-
get_main_info
()¶ Gets the information for the system of interest. Calls
get_state_info()
to obtain information.Returns: Information including the current energy, gradient, and positions for the system of interest Return type: dict
-
get_state_info
(main_info=False, energy=True, positions=True, velocity=False, forces=True, parameters=False, param_deriv=False, periodic_box=False, groups_included=-1)¶ Gets information like the kinetic and potential energy, positions, forces, and topology from an OpenMM state. Some of these may need to be made accessible to user.
Parameters: - simulation (OpenMM simulation object) –
- main_info (bool) – specifies whether to return the topology of the system
- energy (bool) – spcifies whether to get the energy, returned in hartrees(a.u.), default is true.
- positions (bool) – specifies whether to get the positions, returns in nanometers, default is true
- velocity (bool) – specifies whether to get the velocities, default is false
- forces (bool) – specifies whether to get the forces acting on the system, returns as numpy array in jk/mol/nm, as well as the gradients, in au/bohr, default is true
- parameters (bool) – specifies whether to get the parameters of the state, default is false
- param_deriv (bool) – specifies whether to get the parameter derivatives of the state, default is false
- periodic_box (bool) – whether to translate the positions so the center of every molecule lies in the same periodic box, default is false
- groups (list) – a set of indices for which force groups to include when computing forces and energies. Default is all groups
Returns: Information specified by parameters. Keys include ‘energy’, ‘potential’, ‘kinetic’, ‘forces’, ‘gradients’, ‘topology’
Return type: dict
Examples
>>> get_state_info(sim) get_state_info(sim, groups_included=set{0,1,2}) get_state_info(sim, positions=True, forces=True)
-
initialize
(embedding_method)¶ Gets information for the system of interest in its initial state. The simulation object and information dictionary returned by
compute_info()
is saved.Parameters: embedding_method (str) – what QM/MM embedding method to use for initialization. If ‘Mechanical’, all forces are included If ‘Electrostatic’, all coulombic forces are excluded
-
keep_atoms
(atoms)¶ Acts on an OpenMM Modeller object to keep the specified atoms in the MM system and deletes everything else
Parameters: - model (OpenMM Modeller object) –
- atoms (list) – which atoms to keep in an OpenMM Modeller object
Examples
>>> keep_atom(mod, [0,1]) keep_atom(mod, ['O', 'H'])
-
restart
(embedding_method, chkpt_file, restart_forces)¶
-
set_LJ_zero
(OM_system)¶ Removes the Lennard-Jones (van der Waals) force from the system
Parameters: OM_system (OpenMM system object) –
-
set_charge_zero
(OM_system, link_atoms=None)¶ Removes the coulombic forces by setting charges of specified atoms to zero
Parameters: - OM_system (OpenMM system object) –
- link_atoms (list) – link_atoms to set the charge to zero, if link_atoms is None (default), the charge of all particles in the system will be set to zero
Examples
>>> set_charge_zero(system) set_charge_zero(system, link_atoms=[0,1,2])
-
set_external_charges
()¶ Function not implemented for classes
-
set_up_reporters
(simulation)¶ Sets up reporters according to options specified by arguments in md_param. See keywords in the Molecular Dynamics section of the manual for more information.
Parameters: simulation (OpenMM simulation object) – Adds reporters to simulation object
-
take_step
(num)¶ Takes a specified num of steps in the MD simulation
Parameters: num (int) – the number of pure MD steps to take
-
take_updated_step
(force)¶ Updates the system with forces from qmmm and takes a simulation step
Parameters: force (dict) – forces(particle index : forces) in au/bohr to be updated in custom qmmm force and fed into simulation
-
update_forces
(forces, force_obj, simulation)¶ Updates a simulation with external forces
Parameters: - force (dict) – forces(particle index : forces) in au/bohr to be updated in a force object and fed into simulation
- force_obj (OpenMM Force object) – the force object to add the forces to. Can be built in or custom
- simulation (OpenMM simulation object) – where the forces are to be updated in
-
write_pdb
(info)¶ Write a pdb file
Parameters: info (dict) – dictionary containing the topology and positions of the system to write to file