cmp_contact_matrix
- mdtools.structure.cmp_contact_matrix(cm, natms_per_refcmp=1, natms_per_selcmp=1, min_contacts=1)[source]
Convert an
Atom
contact matrix to a compound contact matrix.First sum the contacts of all
Atoms
belonging to the same compound. Then decide based on the number of contacts between theAtoms
of two given compounds whether these two compounds are in contact or not. Two compounds are considered to be in contact if the summed number of contacts between theirAtoms
is equal to or higher than min_contacts.A compound is usually a chemically meaningful subgroup of an
AtomGroup
. This can e.g. be aSegment
,Residue
,fragment
or a singleAtom
. Refer to the MDAnalysis’ user guide for an explanation of these terms. Note that in any case, onlyAtoms
belonging to the originalAtomGroup
are taken into account, even if the compound might comprise additionalAtoms
that are not contained in the originalAtomGroup
.- Parameters:
cm (
array_like
) – (Boolean) contact matrix of shape(m, n)
as e.g. generated bymdtools.structure.contact_matrix()
, wherem
is the number of referenceAtoms
andn
is the number of selectionAtoms
.natms_per_refcmp (
int
orarray_like
, optional) – Number ofAtoms
per reference compound. Can be a single integer or an array of integers. If natms_per_refcmp is a single integer, all reference compounds are assumed to contain the same number ofAtoms
. In this case, natms_per_refcmp must be an integer divisor ofcm.shape[0]
. If natms_per_refcmp is an array of integers, it must contain the number of referenceAtoms
for each single reference compound. In this case,sum(natms_per_refcmp)
must be equal tocm.shape[0]
.natms_per_selcmp (
int
orarray_like
, optional) – Same for selection compounds (natms_per_selcmp is checked againstcm.shape[1]
).min_contacts (
int
, optional) – Two compounds are considered to be in contact if the summed number of contacts between theirAtoms
is equal to or higher than min_contacts. Must be greater than zero. min_contacts is ignored if both natms_per_refcmp and natms_per_selcmp are1
.
- Returns:
cmp_cm (
numpy.ndarray
) – A new (smaller) boolean contact matrix, whose elements evaluate toTrue
if a reference and selection compound have at least min_contacts contacts. If both natms_per_refcmp and natms_per_selcmp are1
, the input contact matrix is returned instead with its data type converted to bool (if it was not already bool before).
See also
mdtools.structure.contact_matrix()
Construct a boolean contact matrix for two MDAnalysis
AtomGroups
mdtools.structure.natms_per_cmp()
Get the number of
Atoms
of each compound in an MDAnalysisAtomGroup
mdtools.structure.cmp_contact_count_matrix()
Take an
Atom
contact matrix and sum the contacts of allAtoms
belonging to the same compound
Notes
Atoms
belonging to the same compound must form a contiguous set in the input contact matrix, otherwise the result will be wrong.If both natms_per_refcmp and natms_per_selcmp are
1
,cmp_contact_matrix(cm)
is equivalent tonp.asarray(cm, dtype=bool)
.Examples
>>> cm = np.array([[ True, True, False, False], ... [ True, False, False, True], ... [False, False, True, True]]) >>> mdt.strc.cmp_contact_matrix(cm, natms_per_refcmp=3) array([[ True, True, True, True]]) >>> mdt.strc.cmp_contact_matrix(cm, natms_per_selcmp=2) array([[ True, False], [ True, True], [False, True]]) >>> mdt.strc.cmp_contact_matrix( ... cm, natms_per_selcmp=2, min_contacts=2 ... ) array([[ True, False], [False, False], [False, True]]) >>> mdt.strc.cmp_contact_matrix( ... cm, natms_per_refcmp=[1, 2], natms_per_selcmp=2 ... ) array([[ True, False], [ True, True]])
If both natms_per_refcmp and natms_per_selcmp are
1
,cmp_contact_matrix(cm)
is equivalent tonp.asarray(cm, dtype=bool)
.>>> mdt.strc.cmp_contact_matrix(cm) is cm True >>> mdt.strc.cmp_contact_matrix(cm, min_contacts=2) is cm True
Edge cases:
>>> cm = np.array([], dtype=bool).reshape(0, 4) >>> mdt.strc.cmp_contact_matrix(cm, natms_per_refcmp=[]) array([], shape=(0, 4), dtype=bool) >>> mdt.strc.cmp_contact_matrix(cm, natms_per_refcmp=1) array([], shape=(0, 4), dtype=bool)
>>> cm = np.array([], dtype=bool).reshape(6, 0) >>> mdt.strc.cmp_contact_matrix( ... cm, natms_per_refcmp=3, natms_per_selcmp=[], min_contacts=2 ... ) array([], shape=(2, 0), dtype=bool)