perses.annihilation.HybridTopologyFactory¶
-
class
perses.annihilation.
HybridTopologyFactory
(topology_proposal, current_positions, new_positions, use_dispersion_correction=False, functions=None, softcore_method='amber', softcore_alpha=None, softcore_beta=None, bond_softening_constant=1.0, angle_softening_constant=1.0, soften_only_new=False)[source]¶ This class generates a hybrid topology based on a perses topology proposal. This class treats atoms in the resulting hybrid system as being from one of four classes:
- unique_old_atom : these atoms are not mapped and only present in the old system. Their interactions will be on for
- lambda=0, off for lambda=1
- unique_new_atom : these atoms are not mapped and only present in the new system. Their interactions will be off
- for lambda=0, on for lambda=1
- core_atom : these atoms are mapped, and are part of a residue that is changing. Their interactions will be those
- corresponding to the old system at lambda=0, and those corresponding to the new system at lambda=1
- environment_atom : these atoms are mapped, and are not part of a changing residue. Their interactions are always
- on and are alchemically unmodified.
Attributes: hybrid_positions
The positions of the hybrid system.
hybrid_system
The hybrid system.
hybrid_topology
An MDTraj hybrid topology for the purpose of writing out trajectories.
new_to_hybrid_atom_map
Give a dictionary that maps new system atoms to the hybrid system.
old_to_hybrid_atom_map
Give a dictionary that maps old system atoms to the hybrid system.
omm_hybrid_topology
An OpenMM format of the hybrid topology.
Methods
handle_harmonic_angles
()This method adds the appropriate interaction for all angles in the hybrid system. handle_harmonic_bonds
()This method adds the appropriate interaction for all bonds in the hybrid system. handle_periodic_torsion_force
()Handle the torsions in the hybrid system in the same way as the angles and bonds. new_positions
(hybrid_positions)Get the positions corresponding to the new system. old_positions
(hybrid_positions)Get the positions corresponding to the old system handle_nonbonded -
__init__
(topology_proposal, current_positions, new_positions, use_dispersion_correction=False, functions=None, softcore_method='amber', softcore_alpha=None, softcore_beta=None, bond_softening_constant=1.0, angle_softening_constant=1.0, soften_only_new=False)[source]¶ Initialize the Hybrid topology factory.
Parameters: - topology_proposal : perses.rjmc.topology_proposal.TopologyProposal object
TopologyProposal object rendered by the ProposalEngine
- current_positions : [n,3] np.ndarray of float
The positions of the “old system”
- new_positions : [m,3] np.ndarray of float
The positions of the “new system”
- use_dispersion_correction : bool, default False
Whether to use the long range correction in the custom sterics force. This is very expensive for NCMC.
- functions : dict, default None
Alchemical functions that determine how each force is scaled with lambda. The keys must be strings with names beginning with lambda_ and ending with each of bonds, angles, torsions, sterics, electrostatics. If functions is none, then the integrator will need to set each of these and parameter derivatives will be unavailable. If functions is not None, all lambdas must be specified.
- softcore_method : str, default ‘default’
The softcore method to use. The options are: default: as an atom is being disappeared, increase softcore strength. For core atoms, don’t use softcore at endpoints, but interpolate to full softcore at lambda=0.5 amber: same as default, but core is excluded from softcore classic: original scheme used by this code. All alchemical atoms interpolate to 0.25 * softcore at lambda=0.5, but don’t use softcore at endpoints.
- softcore_alpha: float, default None
“alpha” parameter of softcore sterics. If None is provided, value will be set to 0.5
- softcore_beta: unit, default None
“beta” parameter of softcore electrostatics. If None is provided, value will be set to 12*unit.angstrom**2
- bond_softening_constant : float
For bonds between unique atoms and unique-core atoms, soften the force constant at the “dummy” endpoint by this factor. If 1.0, do not soften
- angle_softening_constant : float
For bonds between unique atoms and unique-core atoms, soften the force constant at the “dummy” endpoint by this factor. If 1.0, do not soften
Methods
__init__
(topology_proposal, …[, …])Initialize the Hybrid topology factory. handle_harmonic_angles
()This method adds the appropriate interaction for all angles in the hybrid system. handle_harmonic_bonds
()This method adds the appropriate interaction for all bonds in the hybrid system. handle_nonbonded
()handle_periodic_torsion_force
()Handle the torsions in the hybrid system in the same way as the angles and bonds. new_positions
(hybrid_positions)Get the positions corresponding to the new system. old_positions
(hybrid_positions)Get the positions corresponding to the old system Attributes
hybrid_positions
The positions of the hybrid system. hybrid_system
The hybrid system. hybrid_topology
An MDTraj hybrid topology for the purpose of writing out trajectories. new_to_hybrid_atom_map
Give a dictionary that maps new system atoms to the hybrid system. old_to_hybrid_atom_map
Give a dictionary that maps old system atoms to the hybrid system. omm_hybrid_topology
An OpenMM format of the hybrid topology. -
handle_harmonic_angles
()[source]¶ This method adds the appropriate interaction for all angles in the hybrid system. The scheme used, as with bonds, is:
- If the three atoms are all in the core, then we add to the CustomAngleForce and interpolate between the two
- parameters
- Otherwise, we add the angle to a regular angle force.
-
handle_harmonic_bonds
()[source]¶ This method adds the appropriate interaction for all bonds in the hybrid system. The scheme used is:
- If the two atoms are both in the core, then we add to the CustomBondForce and interpolate between the two
- parameters
- Otherwise, we add the bond to a regular bond force.
-
handle_periodic_torsion_force
()[source]¶ Handle the torsions in the hybrid system in the same way as the angles and bonds.
-
hybrid_positions
¶ The positions of the hybrid system. Dimensionality is (n_environment + n_core + n_old_unique + n_new_unique) The positions are assigned by first copying all the mapped positions from the old system in, then copying the mapped positions from the new system.
Returns: - hybrid_positions : [n, 3] Quantity nanometers
-
hybrid_system
¶ The hybrid system.
Returns: - hybrid_system : openmm.System
The system representing a hybrid between old and new topologies
-
hybrid_topology
¶ An MDTraj hybrid topology for the purpose of writing out trajectories. Note that we do not expect this to be able to be parameterized by the openmm forcefield class.
Returns: - hybrid_topology : mdtraj.Topology
-
new_positions
(hybrid_positions)[source]¶ Get the positions corresponding to the new system.
Parameters: - hybrid_positions : [n, 3] np.ndarray with unit
The positions of the hybrid system
Returns: - new_positions : [m, 3] np.ndarray with unit
The positions of the new system
-
new_to_hybrid_atom_map
¶ Give a dictionary that maps new system atoms to the hybrid system.
Returns: - new_to_hybrid_atom_map : dict of {int, int}
The mapping of atoms from the new system to the hybrid
-
old_positions
(hybrid_positions)[source]¶ Get the positions corresponding to the old system
Parameters: - hybrid_positions : [n, 3] np.ndarray with unit
The positions of the hybrid system
Returns: - old_positions : [m, 3] np.ndarray with unit
The positions of the old system
-
old_to_hybrid_atom_map
¶ Give a dictionary that maps old system atoms to the hybrid system.
Returns: - old_to_hybrid_atom_map : dict of {int, int}
The mapping of atoms from the old system to the hybrid
-
omm_hybrid_topology
¶ An OpenMM format of the hybrid topology. Also cannot be used to parameterize system, only to write out trajectories.
Returns: - hybrid_topology : simtk.openmm.app.Topology