API Reference

Network models

class springcraft.GNM(atoms, force_field, masses=None, use_cell_list=True)

Bases: object

This class represents a Gaussian Network Model.

Parameters:
atomsAtomArray, shape=(n,) or ndarray, shape=(n,3), dtype=float

The atoms or their coordinates that are part of the model. It usually contains only CA atoms.

force_fieldForceField, natoms=n

The ForceField that defines the force constants between the given atoms.

massesbool or ndarray, shape=(n,), dtype=float, optional

If an array is given, the Kirchhoff matrix is weighted with the inverse square root of the given masses. If set to true, these masses are automatically inferred from the res_name annotation of atoms, instead. This requires atoms to be an AtomArray. By default no mass-weighting is applied.

use_cell_listbool, optional

If true, a cell list is used to find atoms within cutoff distance instead of checking all pairwise atom distances. This significantly increases the performance for large number of atoms, but is slower for very small systems. If the force_field does not provide a cutoff, no cell list is used regardless.

Attributes:
kirchhoffndarray, shape=(n,n), dtype=float

The Kirchhoff matrix for this model. This is not a copy: Create a copy before modifying this matrix.

covariancendarray, shape=(n*3,n*3), dtype=float

The covariance matrix for this model, i.e. the inverted Kirchhofff matrix. This is not a copy: Create a copy before modifying this matrix.

bfactor(mode_subset=None, tem=None, tem_factors=1.380649e-23)

Computes the isotropic B-factors/temperature factors/ Deby-Waller factors for atoms/coarse-grained beads using the mean-square fluctuation.

These can be used to relate results obtained from ENMs to experimental results.

Parameters:
mode_subsetndarray, shape=(n,), dtype=int, optional

Specifies the subset of modes considered in the MSF computation. Only non-trivial modes can be selected. The first mode is counted as 0 in accordance with Python conventions. If mode_subset is None, all modes except the first trivial mode (0) are included.

temint, float, None, optional

Temperature in Kelvin to compute the temperature scaling factor by multiplying with the Boltzmann constant. If tem is None, no temperature scaling is conducted.

tem_factorsint, float, optional

Factors included in temperature weighting (with K_B as preset).

Returns
——-
bfac_valuesndarray, shape=(n,), dtype=float

B-factors of C-alpha atoms.

dcc(mode_subset=None, norm=True, tem=None, tem_factors=1.380649e-23)

Computes the normalized dynamic cross-correlation between nodes of the GNM.

The DCC is a measure for the correlation in fluctuations exhibited by a given pair of nodes. If normalized, pairs with correlated fluctuations (same phase and period), anticorrelated fluctuations (opposite phase, same period) and non-correlated fluctuations are assigned (normalized) DCC values of 1, -1 and 0 respectively. For results consistent with MSFs, temperature-weighted absolute values can be computed (only relevant if results are not normalized).

Parameters:
mode_subsetndarray, shape=(n,), dtype=int, optional

Specifies the subset of modes considered in the MSF computation. Only non-trivial modes can be selected. The first mode is counted as 0 in accordance with Python conventions. If mode_subset is None, all modes except the first trivial mode (0) are included.

normbool, optional

Normalize the DCC using the MSFs of interacting nodes.

temint, float, None, optional

Temperature in Kelvin to compute the temperature scaling factor by multiplying with the Boltzmann constant. If tem is None, no temperature scaling is conducted.

tem_factorsint, float, optional

Factors included in temperature weighting (with \(k_B\) as preset).

Returns:
dccndarray, shape=(n, n), dtype=float

DCC values for ENM nodes.

Notes

The DCC for a nodepair \(ij\) is computed as:

\[DCC_{ij} = \frac{3 k_B T}{\gamma} \sum_k^L \left[ \frac{\vec{u}_k \cdot \vec{u}_k^T}{\lambda_k} \right]_{ij}\]

with \(\lambda\) and \(\vec{u}\) as Eigenvalues and Eigenvectors corresponding to mode \(k\) of the modeset \(L\).

DCCs can be normalized to MSFs exhibited by two compared nodes following:

\[nDCC_{ij} = \frac{DCC_{ij}}{[DCC_{ii} DCC_{jj}]^{1/2}}\]
eigen()

Compute the Eigenvalues and Eigenvectors of the Kirchhoff matrix.

Returns:
eig_valuesndarray, shape=(k,), dtype=float

Eigenvalues of the Kirchhoff matrix in ascending order.

eig_vectorsndarray, shape=(k,n), dtype=float

Eigenvectors of the Kirchhoff matrix. eig_values[i] corresponds to eigenvectors[i].

frequencies()

Compute the oscillation frequencies of the model.

The first mode corresponds to rigid-body translations/rotations and is omitted in the return value. The returned units are arbitrary and should only be compared relative to each other.

Returns:
frequenciesndarray, shape=(n,), dtype=float

Oscillation frequencies of the model in in ascending order. NaN values mark frequencies corresponding to translations or rotations.

mean_square_fluctuation(mode_subset=None, tem=None, tem_factors=1.380649e-23)

Compute the mean square fluctuation for the atoms according to the GNM. This is equal to the sum of the diagonal of of the GNM covariance matrix, if all k-1 non-trivial modes are considered.

Parameters:
mode_subsetndarray, shape=(n,), dtype=int, optional

Specifies the subset of modes considered in the MSF computation. Only non-trivial modes can be selected. The first mode is counted as 0 in accordance with Python conventions. If mode_subset is None, all modes except the first trivial mode (0) are included.

temint, float, None, optional

Temperature in Kelvin to compute the temperature scaling factor by multiplying with the Boltzmann constant. If tem is None, no temperature scaling is conducted.

tem_factorsint, float, optional

Factors included in temperature weighting (with K_B as preset).

Returns:
msqfndarray, shape=(n,), dtype=float

The mean square fluctuations for each atom in the model.


class springcraft.ANM(atoms, force_field, masses=None, use_cell_list=True)

Bases: object

This class represents an Anisotropic Network Model.

Parameters:
atomsAtomArray, shape=(n,) or ndarray, shape=(n,3), dtype=float

The atoms or their coordinates that are part of the model. It usually contains only CA atoms.

force_fieldForceField, natoms=n

The ForceField that defines the force constants between the given atoms.

massesbool or ndarray, shape=(n,), dtype=float, optional

If an array is given, the Hessian is weighted with the inverse square root of the given masses. If set to true, these masses are automatically inferred from the res_name annotation of atoms, instead. This requires atoms to be an AtomArray. By default no mass-weighting is applied.

use_cell_listbool, optional

If true, a cell list is used to find atoms within cutoff distance instead of checking all pairwise atom distances. This significantly increases the performance for large number of atoms, but is slower for very small systems. If the force_field does not provide a cutoff, no cell list is used regardless.

Attributes:
hessianndarray, shape=(n*3,n*3), dtype=float

The Hessian matrix for this model. Each dimension is partitioned in the form [x1, y1, z1, ... xn, yn, zn]. This is not a copy: Create a copy before modifying this matrix.

covariancendarray, shape=(n*3,n*3), dtype=float

The covariance matrix for this model, i.e. the inverted Hessian. This is not a copy: Create a copy before modifying this matrix.

massesNone or ndarray, shape=(n,), dtype=float

The mass for each atom, None if no mass weighting is applied.

bfactor(mode_subset=None, tem=None, tem_factors=1.380649e-23)

Computes the isotropic B-factors/temperature factors/ Deby-Waller factors for atoms/coarse-grained beads using the mean-square fluctuation. These can be used to relate results obtained from ENMs to experimental results.

Parameters:
mode_subsetndarray, shape=(n,), dtype=int, optional

Specifies the subset of modes considered in the MSF computation. Only non-trivial modes can be selected. The first mode is counted as 0 in accordance with Python conventions. If mode_subset is None, all modes except the first six trivial modes (0-5) are included.

temint, float, None, optional

Temperature in Kelvin to compute the temperature scaling factor by multiplying with the Boltzmann constant. If tem is None, no temperature scaling is conducted.

tem_factorsint, float, optional

Factors included in temperature weighting (with K_B as preset).

Returns
——-
bfac_valuesndarray, shape=(n,), dtype=float

B-factors of C-alpha atoms.

dcc(mode_subset=None, norm=True, tem=None, tem_factors=1.380649e-23)

Computes the normalized dynamic cross-correlation between nodes of the ANM.

The DCC is a measure for the correlation in fluctuations exhibited by a given pair of nodes. If normalized, pairs with correlated fluctuations (same phase and period), anticorrelated fluctuations (opposite phase, same period) and non-correlated fluctuations are assigned (normalized) DCC values of 1, -1 and 0 respectively. For results consistent with MSFs, temperature-weighted absolute values can be computed (only relevant if results are not normalized).

Parameters:
mode_subsetndarray, shape=(n,), dtype=int, optional

Specifies the subset of modes considered in the MSF computation. Only non-trivial modes can be selected. The first mode is counted as 0 in accordance with Python conventions. If mode_subset is None, all modes except the first six trivial modes (0-5) are included.

normbool, optional

Normalize the DCC using the MSFs of interacting nodes.

temint, float, None, optional

Temperature in Kelvin to compute the temperature scaling factor by multiplying with the Boltzmann constant. If tem is None, no temperature scaling is conducted.

tem_factorsint, float, optional

Factors included in temperature weighting (with \(k_B\) as preset).

Returns:
dccndarray, shape=(n, n), dtype=float

DCC values for ENM nodes.

Notes

The DCC for a nodepair \(ij\) is computed as:

\[DCC_{ij} = \frac{3 k_B T}{\gamma} \sum_k^L \left[ \frac{\vec{u}_k \cdot \vec{u}_k^T}{\lambda_k} \right]_{ij}\]

with \(\lambda\) and \(\vec{u}\) as Eigenvalues and Eigenvectors corresponding to mode \(k\) of the modeset \(L\).

DCCs can be normalized to MSFs exhibited by two compared nodes following:

\[nDCC_{ij} = \frac{DCC_{ij}}{[DCC_{ii} DCC_{jj}]^{1/2}}\]
eigen()

Compute the Eigenvalues and Eigenvectors of the Hessian matrix.

The first six Eigenvalues/Eigenvectors correspond to trivial modes (translations/rotations) and are usually omitted in normal mode analysis.

Returns:
eig_valuesndarray, shape=(k,), dtype=float

Eigenvalues of the Hessian matrix in ascending order.

eig_vectorsndarray, shape=(k,n), dtype=float

Eigenvectors of the Hessian matrix. eig_values[i] corresponds to eig_vectors[i].

frequencies()

Computes the frequency associated with each mode.

The first six modes correspond to rigid-body translations/ rotations and are omitted in the return value.

The returned units are arbitrary and should only be compared relative to each other.

Returns:
freqndarray, shape=(n,), dtype=float

The frequency in ascending order of the associated modes’ Eigenvalues.

linear_response(force)

Compute the atom displacement induced by the given force using Linear Response Theory. [1]

Parameters:
forcendarray, shape=(n,3) or shape=(n*3,), dtype=float

The force that is applied to the atoms of the model. The first dimension gives the atom the force is applied on, the second dimension gives the three spatial dimensions. Alternatively, a flattened array in the form [x1, y1, z1, ... xn, yn, zn] can be given.

Returns:
displacementndarray, shape=(n,3), dtype=float

The vector of displacement induced by the given force. The first dimension represents the atom index, the second dimension represents spatial dimension.

References

[1]

M Ikeguchi, J Ueno, M Sato, A Kidera, “Protein Structural Change Upon Ligand Binding: Linear Response Theory.” Phys Rev Lett. 94, 7, 078102 (2005).

mean_square_fluctuation(mode_subset=None, tem=None, tem_factors=1.380649e-23)

Compute the mean square fluctuation for the atoms according to the ANM. This is equal to the sum of the diagonal of each 3x3 superelement of the ANM covariance matrix, if all k-6 non-trivial modes are considered.

Parameters:
mode_subsetndarray, dtype=int, optional

Specifies the subset of modes considered in the MSF computation. Only non-trivial modes can be selected. The first mode is counted as 0 in accordance with Python conventions. If mode_subset is None, all modes except the first six trivial modes (0-5) are included.

temint, float, None, optional

Temperature in Kelvin to compute the temperature scaling factor by multiplying with the Boltzmann constant. If tem is None, no temperature scaling is conducted.

tem_factorsint, float, optional

Factors included in temperature weighting (with K_B as preset).

Returns:
msqfndarray, shape=(n,), dtype=float

The mean square fluctuations for each atom in the model.

normal_mode(index, amplitude, frames, movement='sine')

Create displacements for a trajectory depicting the given normal mode.

This is especially useful for molecular animations of the chosen oscillation mode.

Note, that the first six modes correspond to rigid-body translations/ rotations and are usually omitted in normal mode analysis.

Parameters:
indexint

The index of the oscillation. The index refers to the Eigenvalues obtained from eigen(): Increasing indices refer to oscillations with increasing frequency. The first 6 modes represent tigid body movements (rotations and translations).

amplitudeint

The oscillation amplitude is scaled so that the maximum value for an atom is the given value.

framesint

The number of frames (models) per oscillation.

movement{‘sinusoidal’, ‘triangle’}

Defines how to depict the oscillation. If set to 'sine' the atom movement is sinusoidal. If set to 'triangle' the atom movement is linear with sharp amplitude.

Returns:
displacementndarray, shape=(m,n,3), dtype=float

Atom displacements that depict a single oscillation. m is the number of frames.


Force Fields

class springcraft.ForceField

Bases: object

Subclasses of this abstract base class define the force constants of the modeled springs between atoms in a Elastic network model.

Attributes:
cutoff_distancefloat or None

The interaction of two atoms is only considered, if the distance between them is smaller or equal to this value. If None, the interaction between all atoms is considered.

natomsint or None

The number of atoms in the model. If a ForceField does not depend on the respective atoms, i.e. atom_i and atom_j is unused in force_constant(), this attribute is None instead.

contact_shutdownndarray, shape=(n,), dtype=float, optional

Indices that point to atoms, whose contacts to all other atoms are artificially switched off. If None, no contacts are switched off.

contact_pair_offndarray, shape=(n,2), dtype=int, optional

Indices that point to pairs of atoms, whose contacts are artificially switched off. If None, no contacts are switched off.

contact_pair_onndarray, shape=(n,2), dtype=int, optional

Indices that point to pairs of atoms, whose contacts are are established in any case. If None, no contacts are artificially switched on.

abstract force_constant(atom_i, atom_j, sq_distance)

Get the force constant for the interaction of the given atoms.

ABSTRACT: Override when inheriting.

Parameters:
atom_i, atom_jndarray, shape=(n,), dtype=int

The indices to the first and second atoms in each interacting atom pair.

sq_distancendarray, shape=(n,), dtype=float

The distance between the atoms indicated by atom_i and atom_j.

Notes

Implementations of this method do not need to check whether two atoms are within the cutoff distance of the ForceField: The given pairs of atoms are limited to pairs within cutoff distance of each other. However, if cutoff_distance is None, the atom indices contain the Cartesian product of all atom indices, i.e. each possible combination.


class springcraft.PatchedForceField(force_field, contact_shutdown=None, contact_pair_off=None, contact_pair_on=None, force_constants=None)

Bases: ForceField

This force field wraps another force field and applies custom changes to selected pairs of atoms.

Parameters:
force_fieldForceField

The base force field. For all atoms pairs, that are not patched, the force constant from the base force field is taken

contact_shutdownndarray, shape=(n,), dtype=float, optional

Indices that point to atoms, whose contacts to all other atoms are artificially switched off.

contact_pair_offndarray, shape=(n,2), dtype=int, optional

Indices that point to pairs of atoms, whose contacts are artificially switched off.

contact_pair_onndarray, shape=(n,2), dtype=int, optional

Indices that point to pairs of atoms, whose contacts are are artificially established.

force_constantsndarray, shape=(n,), dtype=float, optional

Individual force constants for artificially established contacts. Must be given, if contact_pair_on is set.


class springcraft.InvariantForceField(cutoff_distance)

Bases: ForceField

This force field treats every interaction with the same force constant.

Parameters:
cutoff_distancefloat

The interaction of two atoms is only considered, if the distance between them is smaller or equal to this value.


class springcraft.HinsenForceField(cutoff_distance=None)

Bases: ForceField

The Hinsen force field was parametrized using the Amber94 force field for a local energy minimum, with crambin as template. In a strict distance-dependent manner, contacts are subdivided into nearest-neighbour pairs along the backbone (r < 4 Å) and mid-/far-range pair interactions (r >= 4 Å). Force constants for these interactions are computed with two distinct formulas. 2.9 Å is the lowest accepted distance between CA atoms. Values below that threshold are set to 2.9 Å.

Parameters:
cutoff_distancefloat, optional

The interaction of two atoms is only considered, if the distance between them is smaller or equal to this value. By default all interactions are included.

References

[1]

K Hinsen et al., “Harmonicity in small proteins.” Chemical Physics 261(1-2): 25-37 (2000).


class springcraft.ParameterFreeForceField(cutoff_distance=None)

Bases: ForceField

The “parameter free ANM” (pfENM) method is an extension of the original implementation of the original ANM forcefield with homogenous parametrization from the Jernigan lab. Unlike in other ANMs, neither distance cutoffs nor distance-dependent spring constants are used. Instead, the residue-pair superelement of the Hessian matrix is weighted by the inverse of the squared distance between residue pairs.

Parameters:
cutoff_distancefloat, optional

The interaction of two atoms is only considered, if the distance between them is smaller or equal to this value. By default all interactions are included.

References

[1]

L Yanga, G Songa, and R L Jernigan, “Protein elastic network models and the ranges of cooperativity.” PNAS. 106, 30, 12347-12352 (2009).


class springcraft.TabulatedForceField(atoms, bonded, intra_chain, inter_chain, cutoff_distance)

Bases: ForceField

This force field uses tabulated force constants for interactions between atoms, based on amino acid type and distance between the atoms.

The distances are separated into bins. A value is within bin[i], if value <= cutoff_distance[i].

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

bonded, intra_chain, inter_chainfloat or ndarray, shape=(k,) or

shape=(20, 20) or shape=(20, 20, k), dtype=float The force constants for interactions between each combination of amino acid type and for each distance bin. The order of amino acids is alphabetically with respect to the one-letter code, i.e. 'ALA', 'CYS', 'ASP', 'GLU', 'PHE', 'GLY', 'HIS', 'ILE', 'LYS', 'LEU', 'MET', 'ASN', 'PRO', 'GLN', 'ARG', 'SER', 'THR', 'VAL', 'TRP' and 'TYR'. bonded gives values for bonded amino acids, intra_chain gives values for non-bonded interactions within the same peptide chain and inter_chain gives values for non-bonded interactions for amino acids in different chains. The possible shapes are:

  • Scalar value: Same value for all amino acid types and distances.

  • 1-dim array: Individual value for each distance bin.

  • 2-dim array: Individual value for each pair of amino acid types. Note the alphabetical order shown above.

  • 3-dim array: Individual value for each distance bin and pair of amino acid types.

cutoff_distancefloat or None or ndarray, shape=(k), dtype=float

If no distance dependent values are given for bonded, intra_chain and inter_chain, this parameter accepts a float, that represents the general cutoff distance, below which interactions between atoms are considered (or None for no cutoff distance). Otherwise, an array of monotonically increasing distance bin edges must be given. The edges represent the right edge of each bin. All interactions at distances above the last edge are not considered.

Attributes:
natomsint or None

The number of atoms in the model.

interaction_matrixndarray, shape=(n, n, k), dtype=float

Force constants between the atoms in atoms. If the tabulated force constants are distance dependent, k is the number of distance bins. Otherwise, k = 1. This is not a copy, modifications on this array affect the force field.

static s_enm_10(atoms)

The sENM10 forcefield by Dehouck and Mikhailov was parametrized by statisctical analysis of a NMR conformational ensemble dataset. Non-bonded interactions between amino acid species are parametrized in an amino acid type-specific manner, with a cutoff distance of 1 nm. Bonded interactions are evaluated with \(10 \, RT/\text{Å}^2\), corresponding to the tenfold mean of all amino acid species interactions at a distance of 3.5 nm.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

Returns:
force_fieldTabulatedForceField

Force field tailored to the sENM10 parameter set.

References

[1]

Y Dehouck, A S Mikhailov, “Effective Harmonic Potentials: Insights into the Internal Cooperativity and Sequence-Specificity of Protein Dynamics.” PLOS Computational Biology 9(8): e1003209 (2013).

static s_enm_13(atoms)

The sENM13 forcefield by Dehouck and Mikhailov was parametrized by statisctical analysis of a NMR conformational ensemble dataset. Non-bonded interactions between amino acid species are parametrized in an amino acid type-specific manner, with a cutoff distance of 1.3 nm. Bonded interactions are evaluated with \(10 \, RT/\text{Å}^2\), corresponding to the tenfold mean of all amino acid species interactions at a distance of 3.5 nm.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

Returns:
force_fieldTabulatedForceField

Force field tailored to the sENM13 parameter set.

References

[1]

Y Dehouck, A S Mikhailov, “Effective Harmonic Potentials: Insights into the Internal Cooperativity and Sequence-Specificity of Protein Dynamics.” PLOS Computational Biology 9(8): e1003209 (2013).

static d_enm(atoms)

The dENM forcefield by Dehouck and Mikhailov was parametrized by statisctical analysis of a NMR conformational ensemble dataset. Non-bonded amino acid interactions are solely assigned depending on the spatial pair distance, ignorant towards interacting amino acid species. Spatial distances are divided into 27 bins. Bonded interactions are evaluated with \(46.83 \, RT/\text{Å}^2\), corresponding to the tenfold mean of all amino acid species interactions at a distance of 3.5 nm.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

Returns:
force_fieldTabulatedForceField

Force field tailored to the dENM parameter set.

References

[1]

Y Dehouck, A S Mikhailov, “Effective Harmonic Potentials: Insights into the Internal Cooperativity and Sequence-Specificity of Protein Dynamics.” PLOS Computational Biology 9(8): e1003209 (2013).

static sd_enm(atoms)

The sdENM forcefield by Dehouck and Mikhailov was parametrized by statistical analysis of a NMR conformational ensemble dataset. Effective harmonic potentials for non-bonded interactions between amino acid pairs are evaluated according to interacting amino acid species as well as the spatial distance between them. Spatial distances are divided into 27 bins, with amino acid specific interaction tables for each distance bin. Bonded interactions are evaluated with \(43.52 \, RT/\text{Å}^2\), corresponding to the tenfold mean of all amino acid species interactions at a distance of 3.5 nm.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

Returns:
force_fieldTabulatedForceField

Force field tailored to the sdENM parameter set.

References

[1]

Y Dehouck, A S Mikhailov, “Effective Harmonic Potentials: Insights into the Internal Cooperativity and Sequence-Specificity of Protein Dynamics.” PLOS Computational Biology 9(8): e1003209 (2013).

static e_anm(atoms, nonbonded_mean=False)

The “extended ANM” (eANM) method discriminates between non-bonded interactions of amino acids within a single polypeptide chain (intrachain) and those present in different chains (interchain) in a residue-specific manner: the former are described by Miyazawa-Jernigan parameters, the latter by Keskin parameters, which are both derived by mean-force statistical analysis of protein structures resolved by X-ray crystallography. Bonded interactions are evaluated with \(82 \, RT/\text{Å}^2\). For noncovalent interactions the cut-off is set to 13 Å.

By averaging over all non-bonded residue-specific parameters, an eANM variant with homogenous parametrization of non-bonded interactions can be derived.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

nonbonded_meanBooleam (optional)

If True, the average of nonbonded interaction tables is computed and used for nonbonded interactions, which yields an homogenous, amino acid-species ignorant parametrization of non-bonded contacts.

Returns:
force_fieldTabulatedForceField

Force field tailored to the eANM method.

References

[1]

K Hamacher, J A McCammon, “Computing the Amino Acid Specificity of Fluctuations in Biomolecular Systems.” J Chem. Theory Comput. 2, 3, 873–878 (2006).

[2]

S Miyazawa, R L Jernigan, “Residue – Residue Potentials with a Favorable Contact Pair Term and an Unfavorable High Packing Density Term, for Simulation and Threading.” J Mol Biol., 256(3) 623-44 (1996).

[3]

O Keskin, I Bahar, R L Jernigan, A Y Badretdinov, O B Ptitsyn, “Empirical solvent-mediated potentials hold for both intra-molecular and inter-molecular inter-residue interactions.” Protein Science, 7 2578-2586 (1998)

e_anm_mj(nonbonded_mean=False)

In this variant of the “extended ANM” (eANM) method, non-bonded interactions between amino acids are parametrized in a residue-specific manner using solely Miyazawa-Jernigan (MJ) parameters. MJ parameters were derived from contact numbers between amino acids in a X-ray structure dataset of globular proteins using the Bethe approximation. Bonded interactions are evaluated with \(82 \, RT/\text{Å}^2\). For noncovalent interactions the cut-off is set to 13 Å.

Averaging over all non-bonded residue-specific parameters yields a homogenous description of non-bonded interactions.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

nonbonded_meanBooleam (optional)

If True, the average of nonbonded interaction tables is computed and used for nonbonded interactions, which yields an homogenous, amino acid-species ignorant parametrization of non-bonded contacts.

Returns:
force_fieldTabulatedForceField

Force field tailored to the eANM method in the MJ parameterset variant.

References

[1]

S Miyazawa, R L Jernigan, “Residue - Residue Potentials with a Favorable Contact Pair Term and an Unfavorable High Packing Density Term, for Simulation and Threading.” J Mol Biol., 256(3) 623-44 (1996).

@TODO Additional citations

e_anm_ke(nonbonded_mean=False)

For this variant of the “extended ANM” (eANM), non-bonded interactions between amino-acid pairs are parametrized in a residue-specific manner using Keskin parameters. This parameterset was derived from contact frequencies between different protein monomers using the methodology established by Miyazawa-Jernigan. Bonded interactions are evaluated with \(82 \, RT/\text{Å}^2\). For noncovalent interactions, the cut-off is set to 13 Å.

Averaging over all non-bonded residue-specific parameters yields a homogenous description of non-bonded interactions.

Parameters:
atomsAtomArray, shape=(n,)

The atoms in the model. Must contain only CA atoms and only canonic amino acids. CA atoms with the same chain ID and adjacent residue IDs are treated as bonded.

nonbonded_meanBooleam (optional)

If True, the average of nonbonded interaction tables is computed and used for nonbonded interactions, which yields an homogenous, amino acid-species ignorant parametrization of non-bonded contacts.

Returns:
force_fieldTabulatedForceField

Force field tailored to the eANM method in the Keskin parameterset variant.

References

[1]

O Keskin, I Bahar, R L Jernigan, A Y Badretdinov, O B Ptitsyn, “Empirical solvent-mediated potentials hold for both intra-molecular and inter-molecular inter-residue interactions.” Protein Science, 7 2578-2586 (1998)


Miscellaneous

springcraft.compute_kirchhoff(coord, force_field, use_cell_list=True)

Compute the Kirchhoff matrix for atoms with given coordinates and the chosen force field.

Parameters:
coordndarray, shape=(n,3), dtype=float

The coordinates.

force_fieldForceField, natoms=n

The ForceField that defines the force constants.

use_cell_listbool, optional

If true, a cell list is used to find atoms within cutoff distance instead of checking all pairwise atom distances. This significantly increases the performance for large number of atoms, but is slower for very small systems. If the force_field does not provide a cutoff, no cell list is used regardless.

Returns:
kirchhoffndarray, shape=(n,n), dtype=float

The computed Kirchhoff matrix.

pairsndarray, shape=(k,2), dtype=int

Indices for interacting atoms, i.e. atoms within cutoff_distance.

springcraft.compute_hessian(coord, force_field, use_cell_list=True)

Compute the Hessian matrix for atoms with given coordinates and the chosen force field.

Parameters:
coordndarray, shape=(n,3), dtype=float

The coordinates.

force_fieldForceField, natoms=n

The ForceField that defines the force constants.

use_cell_listbool, optional

If true, a cell list is used to find atoms within cutoff distance instead of checking all pairwise atom distances. This significantly increases the performance for large number of atoms, but is slower for very small systems. If the force_field does not provide a cutoff, no cell list is used regardless.

Returns:
hessianndarray, shape=(n*3,n*3), dtype=float

The computed Hessian matrix. Each dimension is partitioned in the form [x1, y1, z1, ... xn, yn, zn].

pairsndarray, shape=(k,2), dtype=int

Indices for interacting atoms, i.e. atoms within cutoff_distance.