Geometry API

Classes

AngleMonitor class

class moldesign.AngleMonitor(*atoms)[source]

Bases: moldesign.geom.monitor.Monitor

CONSTRAINT

alias of AngleConstraint

static GETTER(a1, a2, a3)

The angle between bonds a2-a1 and a2-a3

Parameters:a1,a2,a3 (mdt.Atom) – the atoms describing the angle
Returns:the distance
Return type:u.Scalar[length]
static GRAD(a1, a2, a3)

Gradient of the angle between bonds a2-a1 and a2-a3

\[\nabla \theta_{ijkl} = \frac{\partial \theta_{ijkl}}{\partial \mathbf R}\]
Parameters:a1,a2,a3 (mdt.Atom) – the atoms describing the vector

References

https://salilab.org/modeller/9v6/manual/node436.html

NUM_ATOMS = 3
static SETTER(a1, a2, a3, theta, adjustmol=True)

Set the angle between bonds a1-a2 and a3-a2. The atoms will be adjusted along the gradient of the angle. If adjustmol is True and the topology is unambiguous, then the entire molecule’s positions will be modified as well

Parameters:
  • a1,a2,a3 (mdt.Atom) – atoms to adjust
  • theta (u.Scalar[angle]) – new angle to set
  • adjustmol (bool) – Adjust all atoms on either side of this bond?

DihedralMonitor class

class moldesign.DihedralMonitor(*atoms)[source]

Bases: moldesign.geom.monitor.Monitor

CONSTRAINT

alias of DihedralConstraint

static GETTER(a1, a2=None, a3=None, a4=None)

Twist angle of bonds a1-a2 and a4-a3 around around the central bond a2-a3

Can be called as dihedral(a1, a2, a3, a4)
OR dihedral(a2, a2) OR dihedral(bond)
Parameters:
  • a1 (mdt.Bond) – the central bond in the dihedral. OR
  • a1,a2 (mdt.Atom) – the atoms describing the dihedral
  • a3,a4 (mdt.Atom) – (optional) if not passed, a1 and a2 will be treated as the central atoms in this bond, and a3 and a4 will be inferred.
Returns:

angle - [0, 2 pi) radians

Return type:

(units.Scalar[angle])

static GRAD(a1, a2, a3, a4)

Cartesian gradient of a dihedral coordinate,

\[\nabla \theta_{ijkl} = \frac{\partial \theta_{ijkl}}{\partial \mathbf R}\]
Parameters:a1,a2,a3,a4 (mdt.Atom) – the atoms describing the dihedral

References

https://salilab.org/modeller/9v6/manual/node436.html

NUM_ATOMS = 4
static SETTER(a1, a2=None, a3=None, a4=None, theta=None, adjustmol=True)

Set the twist angle of atoms a1 and a4 around the central bond a2-a3. The atoms will be adjusted along the gradient of the angle.

Can be called as set_dihedral(a1, a2, a3, a4, theta, adjustmol=True)
OR set_dihedral(a2, a2, theta, adjustmol=True) OR set_dihedral(bond, theta, adjustmol=True)

If adjustmol is True and the topology is unambiguous, then the entire molecule’s positions will be modified as well

Parameters:
  • a1 (mdt.Bond) – central bond in dihedral
  • a1,a2 (mdt.Atom) – atoms around central bond in dihedral
  • a4 (a3,) –
  • theta (u.Scalar[angle]) – new angle to set
  • adjustmol (bool) – Adjust all atoms on either side of this bond?

DistanceMonitor class

class moldesign.DistanceMonitor(*atoms)[source]

Bases: moldesign.geom.monitor.Monitor

CONSTRAINT

alias of DistanceConstraint

static GETTER(a1, a2)

Return distance between two atoms

Parameters:a1,a2 (mdt.Atom) – the two atoms
Returns:the distance
Return type:u.Scalar[length]
static GRAD(a1, a2)

Gradient of the distance between two atoms,

\[\frac{\partial \mathbf{R}_1}{\partial \mathbf{r}} ||\mathbf{R}_1 - \mathbf{R}_2|| = \frac{\mathbf{R}_1 - \mathbf{R}_2}{||\mathbf{R}_1 - \mathbf{R}_2||}\]
Parameters:a1,a2 (mdt.Atom) – the two atoms
Returns:(gradient w.r.t. first atom, gradient w.r.t. second atom)
Return type:Tuple[u.Vector[length], u.Vector[length]]
NUM_ATOMS = 2
static SETTER(a1, a2, newlength, adjustmol=True)

Set the distance between two atoms. They will be adjusted along the vector separating them. If the two atoms are A) bonded, B) not part of the same ring system, and C) adjustmol is True, then the entire molecule’s positions will be modified as well

Parameters:
  • a1,a2 (mdt.Atom) – atoms to adjust
  • newlength (u.Scalar[length]) – new length to set
  • adjustmol (bool) – Adjust all atoms on either side of this bond?

Functions

angle function

moldesign.angle(a1, a2, a3)[source]

The angle between bonds a2-a1 and a2-a3

Parameters:a1,a2,a3 (mdt.Atom) – the atoms describing the angle
Returns:the distance
Return type:u.Scalar[length]

angle_gradient function

moldesign.angle_gradient(a1, a2, a3)[source]

Gradient of the angle between bonds a2-a1 and a2-a3

\[\nabla \theta_{ijkl} = \frac{\partial \theta_{ijkl}}{\partial \mathbf R}\]
Parameters:a1,a2,a3 (mdt.Atom) – the atoms describing the vector

References

https://salilab.org/modeller/9v6/manual/node436.html

dihedral function

moldesign.dihedral(a1, a2=None, a3=None, a4=None)[source]

Twist angle of bonds a1-a2 and a4-a3 around around the central bond a2-a3

Can be called as dihedral(a1, a2, a3, a4)
OR dihedral(a2, a2) OR dihedral(bond)
Parameters:
  • a1 (mdt.Bond) – the central bond in the dihedral. OR
  • a1,a2 (mdt.Atom) – the atoms describing the dihedral
  • a3,a4 (mdt.Atom) – (optional) if not passed, a1 and a2 will be treated as the central atoms in this bond, and a3 and a4 will be inferred.
Returns:

angle - [0, 2 pi) radians

Return type:

(units.Scalar[angle])

dihedral_gradient function

moldesign.dihedral_gradient(a1, a2, a3, a4)[source]

Cartesian gradient of a dihedral coordinate,

\[\nabla \theta_{ijkl} = \frac{\partial \theta_{ijkl}}{\partial \mathbf R}\]
Parameters:a1,a2,a3,a4 (mdt.Atom) – the atoms describing the dihedral

References

https://salilab.org/modeller/9v6/manual/node436.html

distance function

moldesign.distance(a1, a2)[source]

Return distance between two atoms

Parameters:a1,a2 (mdt.Atom) – the two atoms
Returns:the distance
Return type:u.Scalar[length]

distance_gradient function

moldesign.distance_gradient(a1, a2)[source]

Gradient of the distance between two atoms,

\[\frac{\partial \mathbf{R}_1}{\partial \mathbf{r}} ||\mathbf{R}_1 - \mathbf{R}_2|| = \frac{\mathbf{R}_1 - \mathbf{R}_2}{||\mathbf{R}_1 - \mathbf{R}_2||}\]
Parameters:a1,a2 (mdt.Atom) – the two atoms
Returns:(gradient w.r.t. first atom, gradient w.r.t. second atom)
Return type:Tuple[u.Vector[length], u.Vector[length]]

set_angle function

moldesign.set_angle(a1, a2, a3, theta, adjustmol=True)[source]

Set the angle between bonds a1-a2 and a3-a2. The atoms will be adjusted along the gradient of the angle. If adjustmol is True and the topology is unambiguous, then the entire molecule’s positions will be modified as well

Parameters:
  • a1,a2,a3 (mdt.Atom) – atoms to adjust
  • theta (u.Scalar[angle]) – new angle to set
  • adjustmol (bool) – Adjust all atoms on either side of this bond?

set_dihedral function

moldesign.set_dihedral(a1, a2=None, a3=None, a4=None, theta=None, adjustmol=True)[source]

Set the twist angle of atoms a1 and a4 around the central bond a2-a3. The atoms will be adjusted along the gradient of the angle.

Can be called as set_dihedral(a1, a2, a3, a4, theta, adjustmol=True)
OR set_dihedral(a2, a2, theta, adjustmol=True) OR set_dihedral(bond, theta, adjustmol=True)

If adjustmol is True and the topology is unambiguous, then the entire molecule’s positions will be modified as well

Parameters:
  • a1 (mdt.Bond) – central bond in dihedral
  • a1,a2 (mdt.Atom) – atoms around central bond in dihedral
  • a4 (a3,) –
  • theta (u.Scalar[angle]) – new angle to set
  • adjustmol (bool) – Adjust all atoms on either side of this bond?

set_distance function

moldesign.set_distance(a1, a2, newlength, adjustmol=True)[source]

Set the distance between two atoms. They will be adjusted along the vector separating them. If the two atoms are A) bonded, B) not part of the same ring system, and C) adjustmol is True, then the entire molecule’s positions will be modified as well

Parameters:
  • a1,a2 (mdt.Atom) – atoms to adjust
  • newlength (u.Scalar[length]) – new length to set
  • adjustmol (bool) – Adjust all atoms on either side of this bond?

Unexported API

AngleConstraint class

class moldesign.geom.AngleConstraint(atom1, atom2, atom3, value=None, tolerance=<Quantity(0.0001, 'degree')>, force_constant=<Quantity(1500.0, 'kcalpermol / radian ** 2')>)[source]

Bases: moldesign.geom.constraints.GeometryConstraint

static atomgrad(a1, a2, a3)

Gradient of the angle between bonds a2-a1 and a2-a3

\[\nabla \theta_{ijkl} = \frac{\partial \theta_{ijkl}}{\partial \mathbf R}\]
Parameters:a1,a2,a3 (mdt.Atom) – the atoms describing the vector

References

https://salilab.org/modeller/9v6/manual/node436.html

current()[source]

Return the current value of the constrained quantity

desc = 'angle'
dof = 1
error()[source]
Returns:deviation of current geometry from the constraint
Return type:u.Scalar

DihedralConstraint class

class moldesign.geom.DihedralConstraint(atom1, atom2, atom3, atom4, value=None, tolerance=<Quantity(0.0001, 'degree')>, force_constant=<Quantity(1500.0, 'kcalpermol / radian ** 2')>)[source]

Bases: moldesign.geom.constraints.GeometryConstraint

static atomgrad(a1, a2, a3, a4)

Cartesian gradient of a dihedral coordinate,

\[\nabla \theta_{ijkl} = \frac{\partial \theta_{ijkl}}{\partial \mathbf R}\]
Parameters:a1,a2,a3,a4 (mdt.Atom) – the atoms describing the dihedral

References

https://salilab.org/modeller/9v6/manual/node436.html

current()[source]

Return the current value of the constrained quantity

desc = 'dihedral'
dof = 1
error()[source]
Returns:deviation of current geometry from the constraint
Return type:u.Scalar

DistanceConstraint class

class moldesign.geom.DistanceConstraint(atom1, atom2, value=None, tolerance=<Quantity(1e-06, 'ang')>, force_constant=<Quantity(1000.0, 'kcalpermol / ang ** 2')>)[source]

Bases: moldesign.geom.constraints.GeometryConstraint

static atomgrad(a1, a2)

Gradient of the distance between two atoms,

\[\frac{\partial \mathbf{R}_1}{\partial \mathbf{r}} ||\mathbf{R}_1 - \mathbf{R}_2|| = \frac{\mathbf{R}_1 - \mathbf{R}_2}{||\mathbf{R}_1 - \mathbf{R}_2||}\]
Parameters:a1,a2 (mdt.Atom) – the two atoms
Returns:(gradient w.r.t. first atom, gradient w.r.t. second atom)
Return type:Tuple[u.Vector[length], u.Vector[length]]
current()[source]

Return the current value of the constrained quantity

desc = 'distance'
dof = 1

FixedCoordinate class

class moldesign.geom.FixedCoordinate(atom, vector, value=None, tolerance=<Quantity(1e-06, 'ang')>, force_constant=<Quantity(1000.0, 'kcalpermol / ang ** 2')>)[source]

Bases: moldesign.geom.constraints.GeometryConstraint

Fixes a single, linear degree of freedom for an atom, so that it’s constrained to a plane.

Parameters:
  • atom (moldesign.Atom) – atom to constrain
  • vector (np.ndarray) – direction to constrain
  • value (units.Scalar[length]) – constraint value (i.e., we constrain the dot product of the atom’s position and the normalized direction vector to be this value)
atomgrad(atom=None)[source]
current()[source]

Return the current value of the constrained quantity

desc = 'coordinate'
dof = 1

FixedPosition class

class moldesign.geom.FixedPosition(atom, value=None, tolerance=<Quantity(1e-06, 'ang')>, force_constant=<Quantity(1000.0, 'kcalpermol / ang ** 2')>)[source]

Bases: moldesign.geom.constraints.GeometryConstraint

Fixes a single atom at a given location

Note

The gradient of this function is singular and discontinuous when the constraint is satisfied, leading to poor results in iterative methods such as SHAKE.

In such cases, this constraint should be automatically replaced with three FixedCoordinate constraints on the atom’s x, y, and z coordinates.

atomgrad(atom=None)[source]

Note

For numerical reasons, this returns a vector in the [1,1,1] direction if the constraint is exactly met

Returns:
unit vector from the constraint location to the atom’s actual
location (len=3)
Return type:u.Vector[length]
current()[source]

Return the current value of the constrained quantity

desc = 'position'
dof = 3
error()[source]
Returns:deviation of current geometry from the constraint
Return type:u.Scalar

GeometryConstraint class

class moldesign.geom.GeometryConstraint(atoms, value=None, tolerance=None, force_constant=None)[source]

Bases: object

Base class - Keeps track of a 3D geometry constraint. The constraint is satisfied when self.current() == self.value

Parameters:
  • atoms (List[mdt.Atom]) – atoms involved
  • value (u.Scalar) – constrain the coordinate to this value
  • tolerance (u.Scalar) – absolute tolerance (for iterative constraint enforcement)
  • force_constant (u.Scalar[force]) – optional, only for minimizations and/or use in restraints)
current()[source]

Return the current value of the constrained quantity

desc = 'base constraint'
dof = 1
error()[source]
Returns:deviation of current geometry from the constraint
Return type:u.Scalar
gradient()[source]

Return the gradient of the constrained quantity Requires that self.atomgrad be implemented (or otherwise provided)

\[\nabla G(\mathbf r)\]
restraint_penalty()[source]
Returns:energy penalty for restraints
Return type:u.Scalar[energy]
restraint_penalty_force()[source]
Returns:forces from restraint
Return type:u.Vector[energy]
satisfied()[source]
Returns:True if this constraint is satisfied to within tolerance
Return type:bool

MolecularSymmetry class

class moldesign.geom.MolecularSymmetry(mol, symbol, rms, orientation=None, elems=None, **kwargs)[source]

Bases: object

get_symmetrized_coords(elem)[source]

Symmetrize the molecule based on the symmetry operation This will work as long as the symmetry operation brings each atom closest to a symmetry relation.

Monitor class

class moldesign.geom.Monitor(*atoms)[source]

Bases: object

constrain(**kwargs)[source]

Constrain this coordinate.

This will add a new item to the parent molecule’s constraint list.

Parameters:
  • value (u.Scalar) – constrain the coordinate to this value
  • tolerance (u.Scalar) – absolute tolerance (for iterative constraint enforcement)
  • force_constant (u.Scalar[force]) – optional, only for minimizations and/or use in restraints)
Returns:

the constraint object

Return type:

constraints.GeometryConstraint

gradient()[source]
value

SymmetryElement class

class moldesign.geom.SymmetryElement(symbol, matrix, **kwargs)[source]

Bases: object

Represents an basic, origin-centered symmetry operation: an identity operation (C1/E), inversion center (Ci), mirror plane (Cs), rotation axis (C2,C3,...), or improper rotation (S4, S6, ...)

get_axis()[source]

Returns normal of the plane for Cs or axis of rotation for a Cn :param symm_elem: symmetry element (with attributes ‘symbol’ and ‘matrix’) :return: array of shape (3,)

apply_4x4_transform function

moldesign.geom.apply_4x4_transform(trans, vecs)[source]

Applies a 4x4 transformation vector so one or more 3-D position vector :param trans: :param vecs: :return: transformed position vector

get_symmetry function

moldesign.geom.get_symmetry(mol, tolerance=<Quantity(0.1, 'ang')>, image='symmol', engine=None)

improper_axis_from_matrix function

moldesign.geom.improper_axis_from_matrix(matrix)[source]

Return rotation angle and axis / mirror plane normal from improper rotation matrix.

normalized function

moldesign.geom.normalized(vec)[source]

Return a vector normalized in L2. If vector is 0, return 0

Parameters:vec (u.Vector) – vector to be normalized
Returns:normalized vector
Return type:u.Vector

perpendicular function

moldesign.geom.perpendicular(vec)[source]

safe_arccos function

moldesign.geom.safe_arccos(costheta)[source]

Version of arccos that can handle numerical noise greater than 1.0

shake_positions function

moldesign.geom.shake_positions(mol, prev_positions, max_cycles=100, use_masses=True)[source]

Satisfy all molecular constraints using the SHAKE algorithm

Parameters:
  • mol (moldesign.Molecule) – molecule with unsatisfied constraints
  • prev_positions (mdt.units.Array[length]) – positions prior to move
  • max_cycles (int) – halt and raise an exception after this many cycles

Note

This algorithm does badly if constraint function gradients go to 0.

References

  1. Elber, A. Ruymgaart, B. Hess: SHAKE Parallelization.

    Eur Phys J Spec Top. 2011 Nov 1; 200(1): 211. doi:10.1140/epjst/e2011-01525-9

sub_angles function

moldesign.geom.sub_angles(a, b)[source]

Subtract two angles, keeping the result within [-180,180)