contact_hist_refcmp_diff_selcmp
- mdtools.structure.contact_hist_refcmp_diff_selcmp(cm, natms_per_refcmp=1, natms_per_selcmp=1, minlength=0, dtype=int)[source]
Bin the number of contacts that reference compounds establish to different selection compounds into a histogram.
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
. Alternatively, cm can already be a compound contact matrix as e.g. generated bymdtools.structure.cmp_contact_matrix()
. In this case, you probably want to set natms_per_refcmp and natms_per_selcmp to1
, to keep cm unaltered. It is also possible to parse a compound contact count matrix as e.g. generated bymdtools.structure.cmp_contact_count_matrix()
.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]
).minlength (
int
, optional) – A minimum number of bins for the output array. The output array will have at least this number of elements, though it will be longer if necessary. Seenumpy.bincount()
for further information.dtype (
dtype
, optional) – Data type of the output array.
- Returns:
hist_refcmp_diff_selcmp (
numpy.ndarray
) – Histogram of the number of contacts that reference compounds establish to different selection compounds. Multiple contacts with the same selection compound are not taken into account.
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.contact_hists()
Bin the number of contacts between reference and selection compounds into histograms.
mdtools.structure.contact_hist_refcmp_same_selcmp()
Bin the number of contacts that reference compounds establish to the same selection compound into a histogram.
mdtools.structure.contact_hist_refcmp_selcmp_tot()
Bin the total number of contacts that reference compounds establish to selection compounds into a histogram.
mdtools.structure.contact_hist_refcmp_selcmp_pair()
Bin the number of “bonds” (
Atom
-Atom
contacts) between pairs of reference and selection compoundsnumpy.count_nonzero()
Count the number of non-zero values in an array
numpy.bincount()
Count the number of occurrences of each value in an array of non-negative ints
Notes
Atoms
belonging to the same compound must form a contiguous set in the input contact matrix, otherwise the result will be wrong.About the output array:
- hist_refcmp_diff_selcmp
The first element is the number of reference compounds having no contact with any selection compound, the second element is the number of reference compounds having contact with exactly one selection compound, the third element is the number of reference compounds having contact with exactly two different selection compounds, and so on.
If both natms_per_refcmp and natms_per_selcmp are
1
and cm is a true boolean contact matrix,mdtools.structure.contact_hist_refcmp_diff_selcmp()
andmdtools.structure.contact_hist_refcmp_selcmp_tot()
return the same result.Examples
>>> cm = np.tril(np.ones((5,4), dtype=bool), -1) >>> cm[3, 0] = 0 >>> cm array([[False, False, False, False], [ True, False, False, False], [ True, True, False, False], [False, True, True, False], [ True, True, True, True]]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp(cm) array([1, 1, 2, 0, 1]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp(cm=cm, minlength=7) array([1, 1, 2, 0, 1, 0, 0]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp(cm=cm, minlength=3) array([1, 1, 2, 0, 1]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp(cm=cm, dtype=np.uint32) array([1, 1, 2, 0, 1], dtype=uint32)
>>> mdt.strc.cmp_contact_matrix(cm=cm, natms_per_refcmp=[2, 2, 1]) array([[ True, False, False, False], [ True, True, True, False], [ True, True, True, True]]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp( ... cm=cm, natms_per_refcmp=[2, 2, 1] ... ) array([0, 1, 0, 1, 1])
>>> mdt.strc.cmp_contact_matrix(cm=cm, natms_per_selcmp=2) array([[False, False], [ True, False], [ True, False], [ True, True], [ True, True]]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp( ... cm=cm, natms_per_selcmp=2 ... ) array([1, 2, 2])
>>> mdt.strc.cmp_contact_matrix( ... cm=cm, natms_per_refcmp=[2, 2, 1], natms_per_selcmp=2 ... ) array([[ True, False], [ True, True], [ True, True]]) >>> mdt.strc.contact_hist_refcmp_diff_selcmp( ... cm=cm, natms_per_refcmp=[2, 2, 1], natms_per_selcmp=2 ... ) array([0, 1, 2])
Edge cases:
>>> cm = np.array([], dtype=bool).reshape(0, 4) >>> mdt.strc.contact_hist_refcmp_diff_selcmp( ... cm, natms_per_refcmp=[] ... ) array([], dtype=int64) >>> mdt.strc.contact_hist_refcmp_diff_selcmp(cm, natms_per_refcmp=1) array([], dtype=int64) >>> mdt.strc.contact_hist_refcmp_diff_selcmp( ... cm, natms_per_refcmp=[], minlength=2, dtype=np.uint32 ... ) array([0, 0], dtype=uint32)
>>> cm = np.array([], dtype=bool).reshape(6, 0) >>> mdt.strc.contact_hist_refcmp_diff_selcmp( ... cm, natms_per_refcmp=3, natms_per_selcmp=[] ... ) array([2])